Fork me on GitHub

source: git/classes/DelphesClasses.h@ c749957

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

read xsec from LHE files produced with Wizard

  • Property mode set to 100644
File size: 19.5 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
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
[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
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
[01f9722]77 Float_t CrossSection; // cross-section (read from init, implemented only for Wizard evgen)
[d7d2da3]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, 2)
83};
84
85//---------------------------------------------------------------------------
86
[986d9d5]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
[d7d2da3]98class HepMCEvent: public Event
99{
100public:
101
102 Int_t ProcessID; // unique signal process id | signal_process_id()
[3b465ca]103 Int_t MPI; // number of multi parton interactions | mpi ()
[d7d2da3]104
[59abd43]105 Float_t Weight; // weight for the event
[edeb0f0]106 Float_t CrossSection; // cross-section in pb
107 Float_t CrossSectionError; // cross-section error in pb
[59abd43]108
[d7d2da3]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()
[3b465ca]114 Int_t ID2; // flavour code of second parton | pdf_info()->id2()
[d7d2da3]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
[edeb0f0]124
125
[d7d2da3]126 ClassDef(HepMCEvent, 2)
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
[3b465ca]138
[d7d2da3]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
[151255d]154 Float_t P; // particle momentum
155 Float_t PT; // particle transverse momentum
[d7d2da3]156 Float_t Eta; // particle pseudorapidity
[6d8a29a]157 Float_t Phi; // particle azimuthal angle
158
[d7d2da3]159 Float_t Rapidity; // particle rapidity
[6d8a29a]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
[d7d2da3]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; }
[3b465ca]172
[2b3ef28]173 TLorentzVector P4() const;
[d7d2da3]174
[6d8a29a]175 ClassDef(GenParticle, 2)
[d7d2da3]176};
177
178//---------------------------------------------------------------------------
179
[f59a7b6]180class Vertex: public SortableObject
[d07e957]181{
182public:
[151255d]183
[6d8a29a]184 Float_t T; // vertex position (t component)
[d07e957]185 Float_t X; // vertex position (x component)
186 Float_t Y; // vertex position (y component)
187 Float_t Z; // vertex position (z component)
[5496767]188
[6d8a29a]189 Double_t ErrorT; // vertex position error (t component)
[151255d]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)
[0e2f49b]193
[151255d]194 Int_t Index; // vertex index
195 Int_t NDF; // number of degrees of freedom
[6d8a29a]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
[0e2f49b]203
[5496767]204 TRefArray Constituents; // references to constituents
205
[3c46e17]206 static CompBase *fgCompare; //!
207 const CompBase *GetCompare() const { return fgCompare; }
208
209 ClassDef(Vertex, 3)
[d07e957]210};
211
212//---------------------------------------------------------------------------
213
[d7d2da3]214class MissingET: public TObject
215{
216public:
217 Float_t MET; // mising transverse energy
[4ad7b96]218 Float_t Eta; // mising energy pseudorapidity
[d7d2da3]219 Float_t Phi; // mising energy azimuthal angle
220
[2b3ef28]221 TLorentzVector P4() const;
[4ad7b96]222
[d7d2da3]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
[71648c2]238class Rho: public TObject
239{
240public:
241 Float_t Rho; // rho energy density
[3b465ca]242 Float_t Edges[2]; // pseudorapidity range edges
[71648c2]243
244 ClassDef(Rho, 1)
245};
246
247//---------------------------------------------------------------------------
248
[2e229c9]249class Weight: public TObject
250{
251public:
252 Float_t Weight; // weight for the event
253
254 ClassDef(Weight, 1)
255};
256
257//---------------------------------------------------------------------------
258
[d7d2da3]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
[da00c35]268
[151255d]269 Float_t T; // particle arrival time of flight
[da00c35]270
[d7d2da3]271 Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
272
273 TRefArray Particles; // references to generated particles
274
[151255d]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
[b62c2da]281
[0e0f211]282 Int_t Status; // 1: prompt -- 2: non prompt -- 3: fake
283
[d7d2da3]284 static CompBase *fgCompare; //!
285 const CompBase *GetCompare() const { return fgCompare; }
286
[2b3ef28]287 TLorentzVector P4() const;
[d7d2da3]288
[eb52a5d]289 ClassDef(Photon, 3)
[d7d2da3]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
[da00c35]301
[151255d]302 Float_t T; // particle arrival time of flight
[da00c35]303
[d7d2da3]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
[151255d]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
[b62c2da]316
[d7d2da3]317 static CompBase *fgCompare; //!
318 const CompBase *GetCompare() const { return fgCompare; }
319
[2b3ef28]320 TLorentzVector P4() const;
[d7d2da3]321
[eb52a5d]322 ClassDef(Electron, 3)
[d7d2da3]323};
324
325//---------------------------------------------------------------------------
326
327class Muon: public SortableObject
328{
329public:
330
331 Float_t PT; // muon transverse momentum
332 Float_t Eta; // muon pseudorapidity
333 Float_t Phi; // muon azimuthal angle
334
[151255d]335 Float_t T; // particle arrival time of flight
[da00c35]336
[d7d2da3]337 Int_t Charge; // muon charge
338
339 TRef Particle; // reference to generated particle
340
[151255d]341 Float_t IsolationVar; // isolation variable
342 Float_t IsolationVarRhoCorr; // isolation variable
343 Float_t SumPtCharged; // isolation variable
344 Float_t SumPtNeutral; // isolation variable
345 Float_t SumPtChargedPU; // isolation variable
346 Float_t SumPt; // isolation variable
[b62c2da]347
[d7d2da3]348 static CompBase *fgCompare; //!
349 const CompBase *GetCompare() const { return fgCompare; }
350
[2b3ef28]351 TLorentzVector P4() const;
[d7d2da3]352
[eb52a5d]353 ClassDef(Muon, 3)
[d7d2da3]354};
355
356//---------------------------------------------------------------------------
357
358class Jet: public SortableObject
359{
360public:
361
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
370 Float_t DeltaEta; // jet radius in pseudorapidity
371 Float_t DeltaPhi; // jet radius in azimuthal angle
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
390 Float_t Beta; // (sum pt of charged pile-up constituents)/(sum pt of charged constituents)
391 Float_t BetaStar; // (sum pt of charged constituents coming from hard interaction)/(sum pt of charged constituents)
392 Float_t MeanSqDeltaR; // average distance (squared) between constituent and jet weighted by pt (squared) of constituent
393 Float_t PTD; // average pt between constituent and jet weighted by pt of constituent
394 Float_t FracPt[5]; // (sum pt of constituents within a ring 0.1*i < DeltaR < 0.1*(i+1))/(sum pt of constituents)
[63178fb]395
[666d795]396 Float_t Tau[5]; // N-subjettiness
[edf10ba]397
[ba75867]398 TLorentzVector SoftDroppedJet;
399 TLorentzVector SoftDroppedSubJet1;
400 TLorentzVector SoftDroppedSubJet2;
401
[edf10ba]402 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
403 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
404 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
405
[e9c0d73]406
[edf10ba]407 Int_t NSubJetsTrimmed; // number of subjets trimmed
408 Int_t NSubJetsPruned; // number of subjets pruned
409 Int_t NSubJetsSoftDropped; // number of subjets soft-dropped
410
[e9c0d73]411 Double_t ExclYmerge23;
412 Double_t ExclYmerge34;
413 Double_t ExclYmerge45;
414 Double_t ExclYmerge56;
415
[e4c3fef]416 TRefArray Constituents; // references to constituents
417 TRefArray Particles; // references to generated particles
418
419 static CompBase *fgCompare; //!
420 const CompBase *GetCompare() const { return fgCompare; }
421
[8707eeb]422 TLorentzVector P4() const;
[ba1f1ee]423 TLorentzVector Area;
[24d005f]424
[666d795]425 ClassDef(Jet, 3)
[d7d2da3]426};
427
428//---------------------------------------------------------------------------
429
[3b465ca]430class Track: public SortableObject
[d7d2da3]431{
[3b465ca]432public:
[d7d2da3]433 Int_t PID; // HEP ID number
434
435 Int_t Charge; // track charge
436
[6d8a29a]437 Float_t P; // track momentum
438 Float_t PT; // track transverse momentum
[d7d2da3]439 Float_t Eta; // track pseudorapidity
[6d8a29a]440 Float_t Phi; // track azimuthal angle
441 Float_t CtgTheta; // track cotangent of theta
[5496767]442
[d7d2da3]443 Float_t EtaOuter; // track pseudorapidity at the tracker edge
444 Float_t PhiOuter; // track azimuthal angle at the tracker edge
445
[6d8a29a]446 Float_t T; // track vertex position (t component)
[d7d2da3]447 Float_t X; // track vertex position (x component)
448 Float_t Y; // track vertex position (y component)
449 Float_t Z; // track vertex position (z component)
450
[6d8a29a]451 Float_t TOuter; // track position (t component) at the tracker edge
[d7d2da3]452 Float_t XOuter; // track position (x component) at the tracker edge
453 Float_t YOuter; // track position (y component) at the tracker edge
454 Float_t ZOuter; // track position (z component) at the tracker edge
[e4c3fef]455
[6d8a29a]456 Float_t Xd; // X coordinate of point of closest approach to vertex
457 Float_t Yd; // Y coordinate of point of closest approach to vertex
458 Float_t Zd; // Z coordinate of point of closest approach to vertex
[5496767]459
[6d8a29a]460 Float_t L; // track path length
[151255d]461 Float_t D0; // track transverse impact parameter
462 Float_t DZ; // track longitudinal impact parameter
[5496767]463
[151255d]464 Float_t ErrorP; // track momentum error
[acd0621]465 Float_t ErrorPT; // track transverse momentum error
[151255d]466 Float_t ErrorPhi; // track azimuthal angle error
[6d8a29a]467 Float_t ErrorCtgTheta; // track cotangent of theta error
[5496767]468
[6d8a29a]469 Float_t ErrorT; // time measurement error
470 Float_t ErrorD0; // track transverse impact parameter error
471 Float_t ErrorDZ; // track longitudinal impact parameter error
[d7d2da3]472
473 TRef Particle; // reference to generated particle
474
[2600216]475 Int_t VertexIndex; // reference to vertex
[5496767]476
[d7d2da3]477 static CompBase *fgCompare; //!
478 const CompBase *GetCompare() const { return fgCompare; }
479
[2b3ef28]480 TLorentzVector P4() const;
[d7d2da3]481
[6d8a29a]482 ClassDef(Track, 3)
[d7d2da3]483};
484
485//---------------------------------------------------------------------------
486
[3b465ca]487class Tower: public SortableObject
[d7d2da3]488{
489public:
490 Float_t ET; // calorimeter tower transverse energy
491 Float_t Eta; // calorimeter tower pseudorapidity
492 Float_t Phi; // calorimeter tower azimuthal angle
493
494 Float_t E; // calorimeter tower energy
495
[3db5282]496 Float_t T; // ecal deposit time, averaged by sqrt(EM energy) over all particles, not smeared
[839deb7]497 Int_t NTimeHits; // number of hits contributing to time measurement
[edf10ba]498
[d7d2da3]499 Float_t Eem; // calorimeter tower electromagnetic energy
500 Float_t Ehad; // calorimeter tower hadronic energy
501
502 Float_t Edges[4]; // calorimeter tower edges
503
504 TRefArray Particles; // references to generated particles
505
506 static CompBase *fgCompare; //!
507 const CompBase *GetCompare() const { return fgCompare; }
508
[2b3ef28]509 TLorentzVector P4() const;
[d7d2da3]510
[eb52a5d]511 ClassDef(Tower, 2)
[d7d2da3]512};
513
514//---------------------------------------------------------------------------
515
[8f7db23]516class HectorHit: public SortableObject
517{
518public:
519 Float_t E; // reconstructed energy [GeV]
520
521 Float_t Tx; // angle of the momentum in the horizontal (x,z) plane [urad]
522 Float_t Ty; // angle of the momentum in the verical (y,z) plane [urad]
523
524 Float_t T; // time of flight to the detector [s]
525
526 Float_t X; // horizontal distance to the beam [um]
527 Float_t Y; // vertical distance to the beam [um]
528 Float_t S; // distance to the interaction point [m]
529
[64a4950]530 TRef Particle; // reference to generated particle
531
[8f7db23]532 static CompBase *fgCompare; //!
533 const CompBase *GetCompare() const { return fgCompare; }
534
535 ClassDef(HectorHit, 1)
536};
537
538//---------------------------------------------------------------------------
539
[3b465ca]540class Candidate: public SortableObject
[d7d2da3]541{
542 friend class DelphesFactory;
543
544public:
545 Candidate();
546
547 Int_t PID;
548
549 Int_t Status;
550 Int_t M1, M2, D1, D2;
551
552 Int_t Charge;
553
554 Float_t Mass;
[3b465ca]555
[d7d2da3]556 Int_t IsPU;
[b62c2da]557 Int_t IsRecoPU;
[edf10ba]558
[d7d2da3]559 Int_t IsConstituent;
[5d2481f]560 Int_t IsFromConversion;
[839deb7]561
[fe0273c]562 UInt_t Flavor;
563 UInt_t FlavorAlgo;
564 UInt_t FlavorPhys;
[edf10ba]565
[fe0273c]566 UInt_t BTag;
[edf10ba]567 UInt_t BTagAlgo;
[fe0273c]568 UInt_t BTagPhys;
[edf10ba]569
[264bf40]570 UInt_t TauTag;
[7429c6a]571 Float_t TauWeight;
[d7d2da3]572
573 Float_t Eem;
574 Float_t Ehad;
575
576 Float_t Edges[4];
577 Float_t DeltaEta;
578 Float_t DeltaPhi;
579
[0e2f49b]580 TLorentzVector Momentum, Position, InitialPosition, PositionError, Area;
[80306e6]581
582 Float_t L; // path length
[28c722a]583 Float_t ErrorT; // path length
[80306e6]584 Float_t D0;
585 Float_t ErrorD0;
586 Float_t DZ;
587 Float_t ErrorDZ;
588 Float_t P;
589 Float_t ErrorP;
590 Float_t PT;
591 Float_t ErrorPT;
592 Float_t CtgTheta;
593 Float_t ErrorCtgTheta;
594 Float_t Phi;
595 Float_t ErrorPhi;
[e4c3fef]596
[839deb7]597 Float_t Xd;
598 Float_t Yd;
599 Float_t Zd;
[d7d2da3]600
[a98c7ef]601 // tracking resolution
[5496767]602
[a98c7ef]603 Float_t TrackResolution;
604
[da00c35]605 // PileUpJetID variables
606
[839deb7]607 Int_t NCharged;
608 Int_t NNeutrals;
609 Float_t Beta;
610 Float_t BetaStar;
611 Float_t MeanSqDeltaR;
612 Float_t PTD;
613 Float_t FracPt[5];
[edf10ba]614
[839deb7]615 // Timing information
[edf10ba]616
[839deb7]617 Int_t NTimeHits;
618 std::vector< std::pair< Float_t, Float_t > > ECalEnergyTimePairs;
[e4c3fef]619
[b62c2da]620 // Isolation variables
[edf10ba]621
[b62c2da]622 Float_t IsolationVar;
623 Float_t IsolationVarRhoCorr;
624 Float_t SumPtCharged;
625 Float_t SumPtNeutral;
626 Float_t SumPtChargedPU;
627 Float_t SumPt;
628
[0e2f49b]629 // vertex variables
[5496767]630
[0e2f49b]631 Int_t ClusterIndex;
632 Int_t ClusterNDF;
633 Double_t ClusterSigma;
634 Double_t SumPT2;
635 Double_t BTVSumPT2;
636 Double_t GenDeltaZ;
637 Double_t GenSumPT2;
638
[63178fb]639 // N-subjettiness variables
[e4c3fef]640
641 Float_t Tau[5];
[edf10ba]642
[de6d698]643 // Other Substructure variables
[edf10ba]644
[ba75867]645 TLorentzVector SoftDroppedJet;
646 TLorentzVector SoftDroppedSubJet1;
647 TLorentzVector SoftDroppedSubJet2;
648
[edf10ba]649 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
650 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
651 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
652
653 Int_t NSubJetsTrimmed; // number of subjets trimmed
654 Int_t NSubJetsPruned; // number of subjets pruned
655 Int_t NSubJetsSoftDropped; // number of subjets soft-dropped
[de6d698]656
[e9c0d73]657 // Exclusive clustering variables
658 Double_t ExclYmerge23;
659 Double_t ExclYmerge34;
660 Double_t ExclYmerge45;
661 Double_t ExclYmerge56;
662
[d7d2da3]663 static CompBase *fgCompare; //!
664 const CompBase *GetCompare() const { return fgCompare; }
665
666 void AddCandidate(Candidate *object);
667 TObjArray *GetCandidates();
668
669 Bool_t Overlaps(const Candidate *object) const;
670
671 virtual void Copy(TObject &object) const;
672 virtual TObject *Clone(const char *newname = "") const;
[3b465ca]673 virtual void Clear(Option_t* option = "");
[d7d2da3]674
675private:
676 DelphesFactory *fFactory; //!
677 TObjArray *fArray; //!
[3b465ca]678
[d7d2da3]679 void SetFactory(DelphesFactory *factory) { fFactory = factory; }
680
[3c46e17]681 ClassDef(Candidate, 5)
[d7d2da3]682};
683
684#endif // DelphesClasses_h
685
686
Note: See TracBrowser for help on using the repository browser.