Fork me on GitHub

source: git/classes/DelphesClasses.h@ 4564cba

Last change on this file since 4564cba was cc8716b, checked in by GitHub <noreply@…>, 3 years ago

Update to handle CMS endcap muon detector showers for long-lived particles (#103)

Co-authored-by: christinaw97 <christina.wang@…>

  • Property mode set to 100644
File size: 25.1 KB
RevLine 
[b443089]1/*
2 * Delphes: a framework for fast simulation of a generic collider experiment
3 * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
[1fa50c2]4 *
[b443089]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.
[1fa50c2]9 *
[b443089]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.
[1fa50c2]14 *
[b443089]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
[d7d2da3]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
[341014c]34#include "TLorentzVector.h"
[c18dca6]35#include "TMatrixDSym.h"
[d7d2da3]36#include "TObject.h"
[341014c]37#include "TRef.h"
[d7d2da3]38#include "TRefArray.h"
39
40#include "classes/SortableObject.h"
41
42class DelphesFactory;
43
44//---------------------------------------------------------------------------
45
46class Event: public TObject
47{
48public:
49 Long64_t Number; // event number
50
[151255d]51 Float_t ReadTime; // read time
52 Float_t ProcTime; // processing time
[d7d2da3]53
54 ClassDef(Event, 1)
55};
56
57//---------------------------------------------------------------------------
58
59class LHCOEvent: public Event
60{
61public:
62 Int_t Trigger; // trigger word
63
64 ClassDef(LHCOEvent, 1)
65};
66
67//---------------------------------------------------------------------------
68
69class LHEFEvent: public Event
70{
71public:
72 Int_t ProcessID; // subprocess code for the event | hepup.IDPRUP
73
74 Float_t Weight; // weight for the event | hepup.XWGTUP
[01f9722]75 Float_t CrossSection; // cross-section (read from init, implemented only for Wizard evgen)
[d7d2da3]76 Float_t ScalePDF; // scale in GeV used in the calculation of the PDFs in the event | hepup.SCALUP
77 Float_t AlphaQED; // value of the QED coupling used in the event | hepup.AQEDUP
78 Float_t AlphaQCD; // value of the QCD coupling used in the event | hepup.AQCDUP
79
[769f65b]80 ClassDef(LHEFEvent, 3)
[d7d2da3]81};
82
83//---------------------------------------------------------------------------
84
[986d9d5]85class LHEFWeight: public TObject
86{
87public:
88 Int_t ID; // weight ID
89 Float_t Weight; // weight value
90
91 ClassDef(LHEFWeight, 1)
92};
93
94//---------------------------------------------------------------------------
95
[d7d2da3]96class HepMCEvent: public Event
97{
98public:
99 Int_t ProcessID; // unique signal process id | signal_process_id()
[3b465ca]100 Int_t MPI; // number of multi parton interactions | mpi ()
[d7d2da3]101
[59abd43]102 Float_t Weight; // weight for the event
[edeb0f0]103 Float_t CrossSection; // cross-section in pb
104 Float_t CrossSectionError; // cross-section error in pb
[59abd43]105
[d7d2da3]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()
[3b465ca]111 Int_t ID2; // flavour code of second parton | pdf_info()->id2()
[d7d2da3]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
[769f65b]121 ClassDef(HepMCEvent, 3)
[d7d2da3]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
[3b465ca]133
[d7d2da3]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
[151255d]149 Float_t P; // particle momentum
150 Float_t PT; // particle transverse momentum
[d7d2da3]151 Float_t Eta; // particle pseudorapidity
[6d8a29a]152 Float_t Phi; // particle azimuthal angle
[d7d2da3]153 Float_t Rapidity; // particle rapidity
154
155 Float_t T; // particle vertex position (t component) | hepevt.vhep[number][3]
156 Float_t X; // particle vertex position (x component) | hepevt.vhep[number][0]
157 Float_t Y; // particle vertex position (y component) | hepevt.vhep[number][1]
158 Float_t Z; // particle vertex position (z component) | hepevt.vhep[number][2]
159
[cc8716b]160 Float_t decayX;
161 Float_t decayY;
162 Float_t decayZ;
163 Float_t decayT;
164
[d7d2da3]165 static CompBase *fgCompare; //!
166 const CompBase *GetCompare() const { return fgCompare; }
[3b465ca]167
[2b3ef28]168 TLorentzVector P4() const;
[d7d2da3]169
[6d8a29a]170 ClassDef(GenParticle, 2)
[d7d2da3]171};
172
173//---------------------------------------------------------------------------
174
[f59a7b6]175class Vertex: public SortableObject
[d07e957]176{
177public:
[6d8a29a]178 Float_t T; // vertex position (t component)
[d07e957]179 Float_t X; // vertex position (x component)
180 Float_t Y; // vertex position (y component)
181 Float_t Z; // vertex position (z component)
[5496767]182
[6d8a29a]183 Double_t ErrorT; // vertex position error (t component)
[151255d]184 Double_t ErrorX; // vertex position error (x component)
185 Double_t ErrorY; // vertex position error (y component)
186 Double_t ErrorZ; // vertex position error (z component)
[0e2f49b]187
[151255d]188 Int_t Index; // vertex index
189 Int_t NDF; // number of degrees of freedom
[6d8a29a]190
191 Double_t Sigma; // vertex position (z component) error
192 Double_t SumPT2; // sum pt^2 of tracks attached to the vertex
193 Double_t GenSumPT2; // sum pt^2 of gen tracks attached to the vertex
194
195 Double_t GenDeltaZ; // distance in z to closest generated vertex
196 Double_t BTVSumPT2; // sum pt^2 of tracks attached to the secondary vertex
[0e2f49b]197
[5496767]198 TRefArray Constituents; // references to constituents
199
[3c46e17]200 static CompBase *fgCompare; //!
201 const CompBase *GetCompare() const { return fgCompare; }
202
203 ClassDef(Vertex, 3)
[d07e957]204};
205
206//---------------------------------------------------------------------------
207
[d7d2da3]208class MissingET: public TObject
209{
210public:
211 Float_t MET; // mising transverse energy
[4ad7b96]212 Float_t Eta; // mising energy pseudorapidity
[d7d2da3]213 Float_t Phi; // mising energy azimuthal angle
214
[2b3ef28]215 TLorentzVector P4() const;
[4ad7b96]216
[d7d2da3]217 ClassDef(MissingET, 1)
218};
219
220//---------------------------------------------------------------------------
221
222class ScalarHT: public TObject
223{
224public:
225 Float_t HT; // scalar sum of transverse momenta
226
227 ClassDef(ScalarHT, 1)
228};
229
230//---------------------------------------------------------------------------
231
[71648c2]232class Rho: public TObject
233{
234public:
235 Float_t Rho; // rho energy density
[3b465ca]236 Float_t Edges[2]; // pseudorapidity range edges
[71648c2]237
238 ClassDef(Rho, 1)
239};
240
241//---------------------------------------------------------------------------
242
[2e229c9]243class Weight: public TObject
244{
245public:
246 Float_t Weight; // weight for the event
247
248 ClassDef(Weight, 1)
249};
250
251//---------------------------------------------------------------------------
252
[d7d2da3]253class Photon: public SortableObject
254{
255public:
256 Float_t PT; // photon transverse momentum
257 Float_t Eta; // photon pseudorapidity
258 Float_t Phi; // photon azimuthal angle
259
260 Float_t E; // photon energy
[da00c35]261
[151255d]262 Float_t T; // particle arrival time of flight
[da00c35]263
[d7d2da3]264 Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
265
266 TRefArray Particles; // references to generated particles
267
[151255d]268 Float_t IsolationVar; // isolation variable
269 Float_t IsolationVarRhoCorr; // isolation variable
270 Float_t SumPtCharged; // isolation variable
271 Float_t SumPtNeutral; // isolation variable
272 Float_t SumPtChargedPU; // isolation variable
273 Float_t SumPt; // isolation variable
[b62c2da]274
[0e0f211]275 Int_t Status; // 1: prompt -- 2: non prompt -- 3: fake
276
[d7d2da3]277 static CompBase *fgCompare; //!
278 const CompBase *GetCompare() const { return fgCompare; }
279
[2b3ef28]280 TLorentzVector P4() const;
[d7d2da3]281
[769f65b]282 ClassDef(Photon, 4)
[d7d2da3]283};
284
285//---------------------------------------------------------------------------
286
287class Electron: public SortableObject
288{
289public:
290 Float_t PT; // electron transverse momentum
291 Float_t Eta; // electron pseudorapidity
292 Float_t Phi; // electron azimuthal angle
[da00c35]293
[151255d]294 Float_t T; // particle arrival time of flight
[da00c35]295
[d7d2da3]296 Int_t Charge; // electron charge
297
298 Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
299
300 TRef Particle; // reference to generated particle
301
[151255d]302 Float_t IsolationVar; // isolation variable
303 Float_t IsolationVarRhoCorr; // isolation variable
304 Float_t SumPtCharged; // isolation variable
305 Float_t SumPtNeutral; // isolation variable
306 Float_t SumPtChargedPU; // isolation variable
307 Float_t SumPt; // isolation variable
[b62c2da]308
[0518688]309 Float_t D0; // track transverse impact parameter
310 Float_t DZ; // track longitudinal impact parameter
311 Float_t ErrorD0; // track transverse impact parameter error
312 Float_t ErrorDZ; // track longitudinal impact parameter error
313
[d7d2da3]314 static CompBase *fgCompare; //!
315 const CompBase *GetCompare() const { return fgCompare; }
316
[2b3ef28]317 TLorentzVector P4() const;
[d7d2da3]318
[769f65b]319 ClassDef(Electron, 4)
[d7d2da3]320};
321
322//---------------------------------------------------------------------------
323
324class Muon: public SortableObject
325{
326public:
327 Float_t PT; // muon transverse momentum
328 Float_t Eta; // muon pseudorapidity
329 Float_t Phi; // muon azimuthal angle
330
[151255d]331 Float_t T; // particle arrival time of flight
[da00c35]332
[d7d2da3]333 Int_t Charge; // muon charge
334
335 TRef Particle; // reference to generated particle
336
[151255d]337 Float_t IsolationVar; // isolation variable
338 Float_t IsolationVarRhoCorr; // isolation variable
339 Float_t SumPtCharged; // isolation variable
340 Float_t SumPtNeutral; // isolation variable
341 Float_t SumPtChargedPU; // isolation variable
342 Float_t SumPt; // isolation variable
[b62c2da]343
[0518688]344 Float_t D0; // track transverse impact parameter
345 Float_t DZ; // track longitudinal impact parameter
346 Float_t ErrorD0; // track transverse impact parameter error
347 Float_t ErrorDZ; // track longitudinal impact parameter error
348
[d7d2da3]349 static CompBase *fgCompare; //!
350 const CompBase *GetCompare() const { return fgCompare; }
351
[2b3ef28]352 TLorentzVector P4() const;
[d7d2da3]353
[769f65b]354 ClassDef(Muon, 4)
[d7d2da3]355};
356
357//---------------------------------------------------------------------------
358
359class Jet: public SortableObject
360{
361public:
362 Float_t PT; // jet transverse momentum
363 Float_t Eta; // jet pseudorapidity
364 Float_t Phi; // jet azimuthal angle
365
[da00c35]366 Float_t T; //particle arrival time of flight
367
[d7d2da3]368 Float_t Mass; // jet invariant mass
369
[341014c]370 Float_t DeltaEta; // jet radius in pseudorapidity
371 Float_t DeltaPhi; // jet radius in azimuthal angle
[d7d2da3]372
[151255d]373 UInt_t Flavor; // jet flavor
374 UInt_t FlavorAlgo; // jet flavor
375 UInt_t FlavorPhys; // jet flavor
[edf10ba]376
[fe0273c]377 UInt_t BTag; // 0 or 1 for a jet that has been tagged as containing a heavy quark
[151255d]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
[edf10ba]380
[264bf40]381 UInt_t TauTag; // 0 or 1 for a jet that has been tagged as a tau
[7429c6a]382 Float_t TauWeight; // probability for jet to be identified as tau
[d7d2da3]383
384 Int_t Charge; // tau charge
385
386 Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
387
[edf10ba]388 Int_t NCharged; // number of charged constituents
389 Int_t NNeutrals; // number of neutral constituents
[c614dd7]390
391 Float_t NeutralEnergyFraction; // charged energy fraction
[fd4b326]392 Float_t ChargedEnergyFraction; // neutral energy fraction
[c614dd7]393
[edf10ba]394 Float_t Beta; // (sum pt of charged pile-up constituents)/(sum pt of charged constituents)
395 Float_t BetaStar; // (sum pt of charged constituents coming from hard interaction)/(sum pt of charged constituents)
396 Float_t MeanSqDeltaR; // average distance (squared) between constituent and jet weighted by pt (squared) of constituent
397 Float_t PTD; // average pt between constituent and jet weighted by pt of constituent
398 Float_t FracPt[5]; // (sum pt of constituents within a ring 0.1*i < DeltaR < 0.1*(i+1))/(sum pt of constituents)
[63178fb]399
[666d795]400 Float_t Tau[5]; // N-subjettiness
[edf10ba]401
[ba75867]402 TLorentzVector SoftDroppedJet;
403 TLorentzVector SoftDroppedSubJet1;
404 TLorentzVector SoftDroppedSubJet2;
405
[edf10ba]406 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
407 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
408 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
409
410 Int_t NSubJetsTrimmed; // number of subjets trimmed
411 Int_t NSubJetsPruned; // number of subjets pruned
412 Int_t NSubJetsSoftDropped; // number of subjets soft-dropped
413
[e9c0d73]414 Double_t ExclYmerge23;
415 Double_t ExclYmerge34;
416 Double_t ExclYmerge45;
417 Double_t ExclYmerge56;
[341014c]418
[e4c3fef]419 TRefArray Constituents; // references to constituents
420 TRefArray Particles; // references to generated particles
421
422 static CompBase *fgCompare; //!
423 const CompBase *GetCompare() const { return fgCompare; }
424
[8707eeb]425 TLorentzVector P4() const;
[ba1f1ee]426 TLorentzVector Area;
[24d005f]427
[769f65b]428 ClassDef(Jet, 4)
[d7d2da3]429};
430
431//---------------------------------------------------------------------------
432
[3b465ca]433class Track: public SortableObject
[d7d2da3]434{
[3b465ca]435public:
[d7d2da3]436 Int_t PID; // HEP ID number
437
438 Int_t Charge; // track charge
439
[6d8a29a]440 Float_t P; // track momentum
441 Float_t PT; // track transverse momentum
[d7d2da3]442 Float_t Eta; // track pseudorapidity
[6d8a29a]443 Float_t Phi; // track azimuthal angle
444 Float_t CtgTheta; // track cotangent of theta
[17cd992]445 Float_t C; // track curvature inverse
[fd4b326]446 Float_t Mass; // particle mass
[5496767]447
[d7d2da3]448 Float_t EtaOuter; // track pseudorapidity at the tracker edge
449 Float_t PhiOuter; // track azimuthal angle at the tracker edge
450
[6d8a29a]451 Float_t T; // track vertex position (t component)
[d7d2da3]452 Float_t X; // track vertex position (x component)
453 Float_t Y; // track vertex position (y component)
454 Float_t Z; // track vertex position (z component)
455
[6d8a29a]456 Float_t TOuter; // track position (t component) at the tracker edge
[d7d2da3]457 Float_t XOuter; // track position (x component) at the tracker edge
458 Float_t YOuter; // track position (y component) at the tracker edge
459 Float_t ZOuter; // track position (z component) at the tracker edge
[e4c3fef]460
[6d8a29a]461 Float_t Xd; // X coordinate of point of closest approach to vertex
462 Float_t Yd; // Y coordinate of point of closest approach to vertex
463 Float_t Zd; // Z coordinate of point of closest approach to vertex
[5496767]464
[3a105e5]465 Float_t XFirstHit; // X coordinate of point of closest approach to vertex
466 Float_t YFirstHit; // Y coordinate of point of closest approach to vertex
467 Float_t ZFirstHit; // Z coordinate of point of closest approach to vertex
468
[6d8a29a]469 Float_t L; // track path length
[151255d]470 Float_t D0; // track transverse impact parameter
471 Float_t DZ; // track longitudinal impact parameter
[a95da74]472 Float_t Nclusters; // Number of ionization clusters
[781af69]473 Float_t dNdx; // Number of ionization clusters
[5496767]474
[151255d]475 Float_t ErrorP; // track momentum error
[acd0621]476 Float_t ErrorPT; // track transverse momentum error
[151255d]477 Float_t ErrorPhi; // track azimuthal angle error
[6d8a29a]478 Float_t ErrorCtgTheta; // track cotangent of theta error
[5496767]479
[6d8a29a]480 Float_t ErrorT; // time measurement error
481 Float_t ErrorD0; // track transverse impact parameter error
482 Float_t ErrorDZ; // track longitudinal impact parameter error
[2671df6]483 Float_t ErrorC; // track curvature error
484
485 // track covariance off-diagonal terms
[fd4b326]486 Float_t ErrorD0Phi;
487 Float_t ErrorD0C;
488 Float_t ErrorD0DZ;
489 Float_t ErrorD0CtgTheta;
490 Float_t ErrorPhiC;
491 Float_t ErrorPhiDZ;
492 Float_t ErrorPhiCtgTheta ;
493 Float_t ErrorCDZ;
494 Float_t ErrorCCtgTheta;
495 Float_t ErrorDZCtgTheta;
[d7d2da3]496
497 TRef Particle; // reference to generated particle
498
[2600216]499 Int_t VertexIndex; // reference to vertex
[5496767]500
[d7d2da3]501 static CompBase *fgCompare; //!
502 const CompBase *GetCompare() const { return fgCompare; }
503
[2b3ef28]504 TLorentzVector P4() const;
[2671df6]505 TMatrixDSym CovarianceMatrix() const;
[d7d2da3]506
[6d8a29a]507 ClassDef(Track, 3)
[d7d2da3]508};
509
510//---------------------------------------------------------------------------
511
[3b465ca]512class Tower: public SortableObject
[d7d2da3]513{
514public:
515 Float_t ET; // calorimeter tower transverse energy
516 Float_t Eta; // calorimeter tower pseudorapidity
517 Float_t Phi; // calorimeter tower azimuthal angle
518
519 Float_t E; // calorimeter tower energy
520
[3db5282]521 Float_t T; // ecal deposit time, averaged by sqrt(EM energy) over all particles, not smeared
[839deb7]522 Int_t NTimeHits; // number of hits contributing to time measurement
[edf10ba]523
[d7d2da3]524 Float_t Eem; // calorimeter tower electromagnetic energy
525 Float_t Ehad; // calorimeter tower hadronic energy
[61dccd3]526 Float_t Etrk; // total charged energy hitting tower
[d7d2da3]527
528 Float_t Edges[4]; // calorimeter tower edges
529
530 TRefArray Particles; // references to generated particles
531
532 static CompBase *fgCompare; //!
533 const CompBase *GetCompare() const { return fgCompare; }
534
[2b3ef28]535 TLorentzVector P4() const;
[d7d2da3]536
[61dccd3]537 ClassDef(Tower, 3)
[d7d2da3]538};
539
540//---------------------------------------------------------------------------
541
[ededa33]542class ParticleFlowCandidate: public SortableObject
543{
544
545public:
546 Int_t PID; // HEP ID number
547
548 Int_t Charge; // track charge
549
550 Float_t E; // reconstructed energy [GeV]
551 Float_t P; // track momentum
552 Float_t PT; // track transverse momentum
553 Float_t Eta; // track pseudorapidity
554 Float_t Phi; // track azimuthal angle
555 Float_t CtgTheta; // track cotangent of theta
[17cd992]556 Float_t C; // track curvature inverse
[fd4b326]557 Float_t Mass; // particle mass
[ededa33]558
559 Float_t EtaOuter; // track pseudorapidity at the tracker edge
560 Float_t PhiOuter; // track azimuthal angle at the tracker edge
561
562 Float_t T; // track vertex position (t component)
563 Float_t X; // track vertex position (x component)
564 Float_t Y; // track vertex position (y component)
565 Float_t Z; // track vertex position (z component)
566
567 Float_t TOuter; // track position (t component) at the tracker edge
568 Float_t XOuter; // track position (x component) at the tracker edge
569 Float_t YOuter; // track position (y component) at the tracker edge
570 Float_t ZOuter; // track position (z component) at the tracker edge
571
572 Float_t Xd; // X coordinate of point of closest approach to vertex
573 Float_t Yd; // Y coordinate of point of closest approach to vertex
574 Float_t Zd; // Z coordinate of point of closest approach to vertex
575
[3a105e5]576 Float_t XFirstHit; // X coordinate of point of closest approach to vertex
577 Float_t YFirstHit; // Y coordinate of point of closest approach to vertex
578 Float_t ZFirstHit; // Z coordinate of point of closest approach to vertex
579
[ededa33]580 Float_t L; // track path length
581 Float_t D0; // track transverse impact parameter
582 Float_t DZ; // track longitudinal impact parameter
[a95da74]583 Float_t Nclusters; // Number of ionization clusters
[781af69]584 Float_t dNdx; // Number of ionization clusters
[ededa33]585
586 Float_t ErrorP; // track momentum error
587 Float_t ErrorPT; // track transverse momentum error
588 Float_t ErrorPhi; // track azimuthal angle error
589 Float_t ErrorCtgTheta; // track cotangent of theta error
590
591 Float_t ErrorT; // time measurement error
592 Float_t ErrorD0; // track transverse impact parameter error
593 Float_t ErrorDZ; // track longitudinal impact parameter error
[2671df6]594 Float_t ErrorC; // track curvature error
595
596 // track covariance off-diagonal terms
[fd4b326]597 Float_t ErrorD0Phi;
598 Float_t ErrorD0C;
599 Float_t ErrorD0DZ;
600 Float_t ErrorD0CtgTheta;
601 Float_t ErrorPhiC;
602 Float_t ErrorPhiDZ;
603 Float_t ErrorPhiCtgTheta ;
604 Float_t ErrorCDZ;
605 Float_t ErrorCCtgTheta;
606 Float_t ErrorDZCtgTheta;
[ededa33]607
608 Int_t VertexIndex; // reference to vertex
609
610 static CompBase *fgCompare; //!
611 const CompBase *GetCompare() const { return fgCompare; }
612
613 TLorentzVector P4() const;
[2671df6]614 TMatrixDSym CovarianceMatrix() const;
[ededa33]615
616 Int_t NTimeHits; // number of hits contributing to time measurement
617
618 Float_t Eem; // calorimeter tower electromagnetic energy
619 Float_t Ehad; // calorimeter tower hadronic energy
[61dccd3]620 Float_t Etrk; // total charged energy hitting tower
[ededa33]621
622 Float_t Edges[4]; // calorimeter tower edges
623
624 TRefArray Particles; // references to generated particles
625
[61dccd3]626 ClassDef(ParticleFlowCandidate, 3)
[ededa33]627
628};
629
630//---------------------------------------------------------------------------
631
[8f7db23]632class HectorHit: public SortableObject
633{
634public:
635 Float_t E; // reconstructed energy [GeV]
636
637 Float_t Tx; // angle of the momentum in the horizontal (x,z) plane [urad]
638 Float_t Ty; // angle of the momentum in the verical (y,z) plane [urad]
639
640 Float_t T; // time of flight to the detector [s]
641
642 Float_t X; // horizontal distance to the beam [um]
643 Float_t Y; // vertical distance to the beam [um]
644 Float_t S; // distance to the interaction point [m]
645
[64a4950]646 TRef Particle; // reference to generated particle
647
[8f7db23]648 static CompBase *fgCompare; //!
649 const CompBase *GetCompare() const { return fgCompare; }
650
651 ClassDef(HectorHit, 1)
652};
[cc8716b]653//---------------------------------------------------------------------------
654
655class CscCluster: public SortableObject
656{
657public:
658 Float_t Eta; // eta of LLP
659 Float_t Phi; // phi of LLP
660 Float_t PT; // pt of LLP
661 Float_t Px;// px of LLP
662 Float_t Py;// py of LLP
663 Float_t Pz;// pz of LLP
664 Float_t E; // E of LLP
665 Float_t Ehad; // had energy of LLP
666 Float_t Eem; // em energy of LLP
667 Float_t pid; // LLP pid
668 Float_t T; // LLP decay time-photon travel time
669 Float_t X; // LLP decay x
670 Float_t Y; // LLP decay y
671 Float_t Z; // LLP decay z
672 Float_t R; // LLP decay z
673 Float_t beta; // LLP beta
674 Float_t ctau; //LLP ctau
675
676
677 static CompBase *fgCompare; //!
678 const CompBase *GetCompare() const { return fgCompare; }
679
680 ClassDef(CscCluster, 5)
681};
[8f7db23]682
683//---------------------------------------------------------------------------
684
[3b465ca]685class Candidate: public SortableObject
[d7d2da3]686{
687 friend class DelphesFactory;
688
689public:
690 Candidate();
691
692 Int_t PID;
693
694 Int_t Status;
695 Int_t M1, M2, D1, D2;
696
697 Int_t Charge;
698
699 Float_t Mass;
[3b465ca]700
[d7d2da3]701 Int_t IsPU;
[b62c2da]702 Int_t IsRecoPU;
[edf10ba]703
[d7d2da3]704 Int_t IsConstituent;
[5d2481f]705 Int_t IsFromConversion;
[839deb7]706
[fe0273c]707 UInt_t Flavor;
708 UInt_t FlavorAlgo;
709 UInt_t FlavorPhys;
[edf10ba]710
[fe0273c]711 UInt_t BTag;
[edf10ba]712 UInt_t BTagAlgo;
[fe0273c]713 UInt_t BTagPhys;
[edf10ba]714
[264bf40]715 UInt_t TauTag;
[7429c6a]716 Float_t TauWeight;
[d7d2da3]717
718 Float_t Eem;
719 Float_t Ehad;
[61dccd3]720 Float_t Etrk;
[d7d2da3]721
722 Float_t Edges[4];
723 Float_t DeltaEta;
724 Float_t DeltaPhi;
725
[cc8716b]726 TLorentzVector Momentum, Position, InitialPosition, PositionError, DecayPosition, Area;
[80306e6]727
728 Float_t L; // path length
[3051ea17]729 Float_t DZ;
730 Float_t ErrorDZ;
[28c722a]731 Float_t ErrorT; // path length
[80306e6]732 Float_t D0;
733 Float_t ErrorD0;
[3051ea17]734 Float_t C;
735 Float_t ErrorC;
[80306e6]736 Float_t P;
737 Float_t ErrorP;
738 Float_t PT;
739 Float_t ErrorPT;
740 Float_t CtgTheta;
741 Float_t ErrorCtgTheta;
742 Float_t Phi;
743 Float_t ErrorPhi;
[e4c3fef]744
[a95da74]745 Float_t Nclusters; // Number of ionization clusters
[781af69]746 Float_t dNdx; // Number of ionization clusters per unit length
[a95da74]747
[839deb7]748 Float_t Xd;
749 Float_t Yd;
750 Float_t Zd;
[d7d2da3]751
[3a105e5]752 Float_t XFirstHit;
753 Float_t YFirstHit;
754 Float_t ZFirstHit;
755
[a98c7ef]756 // tracking resolution
[5496767]757
[a98c7ef]758 Float_t TrackResolution;
759
[da00c35]760 // PileUpJetID variables
761
[839deb7]762 Int_t NCharged;
763 Int_t NNeutrals;
764 Float_t Beta;
765 Float_t BetaStar;
766 Float_t MeanSqDeltaR;
767 Float_t PTD;
768 Float_t FracPt[5];
[c614dd7]769 Float_t NeutralEnergyFraction; // charged energy fraction
[fd4b326]770 Float_t ChargedEnergyFraction; // neutral energy fraction
[c614dd7]771
[edf10ba]772
[839deb7]773 // Timing information
[edf10ba]774
[839deb7]775 Int_t NTimeHits;
[77e9ae1]776 std::vector<std::pair<Float_t, Float_t> > ECalEnergyTimePairs;
[e4c3fef]777
[b62c2da]778 // Isolation variables
[edf10ba]779
[b62c2da]780 Float_t IsolationVar;
781 Float_t IsolationVarRhoCorr;
782 Float_t SumPtCharged;
783 Float_t SumPtNeutral;
784 Float_t SumPtChargedPU;
785 Float_t SumPt;
786
[3051ea17]787 // ACTS compliant 6x6 track covariance (D0, phi, Curvature, dz, ctg(theta))
[c18dca6]788
[3051ea17]789 TMatrixDSym TrackCovariance;
[c18dca6]790
[0e2f49b]791 // vertex variables
[5496767]792
[0e2f49b]793 Int_t ClusterIndex;
794 Int_t ClusterNDF;
795 Double_t ClusterSigma;
796 Double_t SumPT2;
797 Double_t BTVSumPT2;
798 Double_t GenDeltaZ;
799 Double_t GenSumPT2;
800
[63178fb]801 // N-subjettiness variables
[e4c3fef]802
803 Float_t Tau[5];
[edf10ba]804
[de6d698]805 // Other Substructure variables
[edf10ba]806
[ba75867]807 TLorentzVector SoftDroppedJet;
808 TLorentzVector SoftDroppedSubJet1;
809 TLorentzVector SoftDroppedSubJet2;
810
[edf10ba]811 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
812 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
813 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
814
815 Int_t NSubJetsTrimmed; // number of subjets trimmed
816 Int_t NSubJetsPruned; // number of subjets pruned
817 Int_t NSubJetsSoftDropped; // number of subjets soft-dropped
[de6d698]818
[e9c0d73]819 // Exclusive clustering variables
820 Double_t ExclYmerge23;
821 Double_t ExclYmerge34;
822 Double_t ExclYmerge45;
823 Double_t ExclYmerge56;
[341014c]824
[7e83689]825 // event characteristics variables
826 Double_t ParticleDensity; // particle multiplicity density in the proximity of the particle
[fd4b326]827
[d7d2da3]828 static CompBase *fgCompare; //!
829 const CompBase *GetCompare() const { return fgCompare; }
830
831 void AddCandidate(Candidate *object);
832 TObjArray *GetCandidates();
833
834 Bool_t Overlaps(const Candidate *object) const;
835
836 virtual void Copy(TObject &object) const;
837 virtual TObject *Clone(const char *newname = "") const;
[341014c]838 virtual void Clear(Option_t *option = "");
[d7d2da3]839
840private:
841 DelphesFactory *fFactory; //!
842 TObjArray *fArray; //!
[3b465ca]843
[d7d2da3]844 void SetFactory(DelphesFactory *factory) { fFactory = factory; }
845
[769f65b]846 ClassDef(Candidate, 6)
[d7d2da3]847};
848
849#endif // DelphesClasses_h
Note: See TracBrowser for help on using the repository browser.