Fork me on GitHub

source: git/classes/DelphesClasses.h@ af6ce6a

Timing
Last change on this file since af6ce6a was 84a097e, checked in by Michele Selvaggi <michele.selvaggi@…>, 5 years ago

first iteration of adding timing in Delphes

  • Property mode set to 100644
File size: 22.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 "TLorentzVector.h"
35#include "TObject.h"
36#include "TRef.h"
37#include "TRefArray.h"
38
39#include "classes/SortableObject.h"
40
41class DelphesFactory;
42
43//---------------------------------------------------------------------------
44
45class Event: public TObject
46{
47public:
48 Long64_t Number; // event number
49
50 Float_t ReadTime; // read time
51 Float_t ProcTime; // processing time
52
53 ClassDef(Event, 1)
54};
55
56//---------------------------------------------------------------------------
57
58class LHCOEvent: public Event
59{
60public:
61 Int_t Trigger; // trigger word
62
63 ClassDef(LHCOEvent, 1)
64};
65
66//---------------------------------------------------------------------------
67
68class LHEFEvent: public Event
69{
70public:
71 Int_t ProcessID; // subprocess code for the event | hepup.IDPRUP
72
73 Float_t Weight; // weight for the event | hepup.XWGTUP
74 Float_t CrossSection; // cross-section (read from init, implemented only for Wizard evgen)
75 Float_t ScalePDF; // scale in GeV used in the calculation of the PDFs in the event | hepup.SCALUP
76 Float_t AlphaQED; // value of the QED coupling used in the event | hepup.AQEDUP
77 Float_t AlphaQCD; // value of the QCD coupling used in the event | hepup.AQCDUP
78
79 ClassDef(LHEFEvent, 3)
80};
81
82//---------------------------------------------------------------------------
83
84class LHEFWeight: public TObject
85{
86public:
87 Int_t ID; // weight ID
88 Float_t Weight; // weight value
89
90 ClassDef(LHEFWeight, 1)
91};
92
93//---------------------------------------------------------------------------
94
95class HepMCEvent: public Event
96{
97public:
98 Int_t ProcessID; // unique signal process id | signal_process_id()
99 Int_t MPI; // number of multi parton interactions | mpi ()
100
101 Float_t Weight; // weight for the event
102 Float_t CrossSection; // cross-section in pb
103 Float_t CrossSectionError; // cross-section error in pb
104
105 Float_t Scale; // energy scale, see hep-ph/0109068 | event_scale()
106 Float_t AlphaQED; // QED coupling, see hep-ph/0109068 | alphaQED()
107 Float_t AlphaQCD; // QCD coupling, see hep-ph/0109068 | alphaQCD()
108
109 Int_t ID1; // flavour code of first parton | pdf_info()->id1()
110 Int_t ID2; // flavour code of second parton | pdf_info()->id2()
111
112 Float_t X1; // fraction of beam momentum carried by first parton ("beam side") | pdf_info()->x1()
113 Float_t X2; // fraction of beam momentum carried by second parton ("target side") | pdf_info()->x2()
114
115 Float_t ScalePDF; // Q-scale used in evaluation of PDF's (in GeV) | pdf_info()->scalePDF()
116
117 Float_t PDF1; // PDF (id1, x1, Q) | pdf_info()->pdf1()
118 Float_t PDF2; // PDF (id2, x2, Q) | pdf_info()->pdf2()
119
120 ClassDef(HepMCEvent, 3)
121};
122
123//---------------------------------------------------------------------------
124
125class GenParticle: public SortableObject
126{
127public:
128 Int_t PID; // particle HEP ID number | hepevt.idhep[number]
129
130 Int_t Status; // particle status | hepevt.isthep[number]
131 Int_t IsPU; // 0 or 1 for particles from pile-up interactions
132
133 Int_t M1; // particle 1st mother | hepevt.jmohep[number][0] - 1
134 Int_t M2; // particle 2nd mother | hepevt.jmohep[number][1] - 1
135
136 Int_t D1; // particle 1st daughter | hepevt.jdahep[number][0] - 1
137 Int_t D2; // particle last daughter | hepevt.jdahep[number][1] - 1
138
139 Int_t Charge; // particle charge
140
141 Float_t Mass; // particle mass
142
143 Float_t E; // particle energy | hepevt.phep[number][3]
144 Float_t Px; // particle momentum vector (x component) | hepevt.phep[number][0]
145 Float_t Py; // particle momentum vector (y component) | hepevt.phep[number][1]
146 Float_t Pz; // particle momentum vector (z component) | hepevt.phep[number][2]
147
148 Float_t P; // particle momentum
149 Float_t PT; // particle transverse momentum
150 Float_t Eta; // particle pseudorapidity
151 Float_t Phi; // particle azimuthal angle
152
153 Float_t Rapidity; // particle rapidity
154 Float_t CtgTheta; // particle cotangent of theta
155
156 Float_t D0; // particle transverse impact parameter
157 Float_t DZ; // particle longitudinal impact parameter
158
159 Float_t T; // particle vertex position (t component) | hepevt.vhep[number][3]
160 Float_t X; // particle vertex position (x component) | hepevt.vhep[number][0]
161 Float_t Y; // particle vertex position (y component) | hepevt.vhep[number][1]
162 Float_t Z; // particle vertex position (z component) | hepevt.vhep[number][2]
163
164 static CompBase *fgCompare; //!
165 const CompBase *GetCompare() const { return fgCompare; }
166
167 TLorentzVector P4() const;
168
169 ClassDef(GenParticle, 2)
170};
171
172//---------------------------------------------------------------------------
173
174class Vertex: public SortableObject
175{
176public:
177 Float_t T; // vertex position (t component)
178 Float_t X; // vertex position (x component)
179 Float_t Y; // vertex position (y component)
180 Float_t Z; // vertex position (z component)
181
182 Double_t ErrorT; // vertex position error (t component)
183 Double_t ErrorX; // vertex position error (x component)
184 Double_t ErrorY; // vertex position error (y component)
185 Double_t ErrorZ; // vertex position error (z component)
186
187 Int_t Index; // vertex index
188 Int_t NDF; // number of degrees of freedom
189
190 Double_t Sigma; // vertex position (z component) error
191 Double_t SumPT2; // sum pt^2 of tracks attached to the vertex
192 Double_t GenSumPT2; // sum pt^2 of gen tracks attached to the vertex
193
194 Double_t GenDeltaZ; // distance in z to closest generated vertex
195 Double_t BTVSumPT2; // sum pt^2 of tracks attached to the secondary vertex
196
197 TRefArray Constituents; // references to constituents
198
199 static CompBase *fgCompare; //!
200 const CompBase *GetCompare() const { return fgCompare; }
201
202 ClassDef(Vertex, 3)
203};
204
205//---------------------------------------------------------------------------
206
207class MissingET: public TObject
208{
209public:
210 Float_t MET; // mising transverse energy
211 Float_t Eta; // mising energy pseudorapidity
212 Float_t Phi; // mising energy azimuthal angle
213
214 TLorentzVector P4() const;
215
216 ClassDef(MissingET, 1)
217};
218
219//---------------------------------------------------------------------------
220
221class ScalarHT: public TObject
222{
223public:
224 Float_t HT; // scalar sum of transverse momenta
225
226 ClassDef(ScalarHT, 1)
227};
228
229//---------------------------------------------------------------------------
230
231class Rho: public TObject
232{
233public:
234 Float_t Rho; // rho energy density
235 Float_t Edges[2]; // pseudorapidity range edges
236
237 ClassDef(Rho, 1)
238};
239
240//---------------------------------------------------------------------------
241
242class Weight: public TObject
243{
244public:
245 Float_t Weight; // weight for the event
246
247 ClassDef(Weight, 1)
248};
249
250//---------------------------------------------------------------------------
251
252class Photon: public SortableObject
253{
254public:
255 Float_t PT; // photon transverse momentum
256 Float_t Eta; // photon pseudorapidity
257 Float_t Phi; // photon azimuthal angle
258
259 Float_t E; // photon energy
260
261 Float_t T; // particle arrival time of flight
262
263 Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
264
265 TRefArray Particles; // references to generated particles
266
267 Float_t IsolationVar; // isolation variable
268 Float_t IsolationVarRhoCorr; // isolation variable
269 Float_t SumPtCharged; // isolation variable
270 Float_t SumPtNeutral; // isolation variable
271 Float_t SumPtChargedPU; // isolation variable
272 Float_t SumPt; // isolation variable
273
274 Int_t Status; // 1: prompt -- 2: non prompt -- 3: fake
275
276 static CompBase *fgCompare; //!
277 const CompBase *GetCompare() const { return fgCompare; }
278
279 TLorentzVector P4() const;
280
281 ClassDef(Photon, 4)
282};
283
284//---------------------------------------------------------------------------
285
286class Electron: public SortableObject
287{
288public:
289 Float_t PT; // electron transverse momentum
290 Float_t Eta; // electron pseudorapidity
291 Float_t Phi; // electron azimuthal angle
292
293 Float_t T; // particle arrival time of flight
294
295 Int_t Charge; // electron charge
296
297 Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
298
299 TRef Particle; // reference to generated particle
300
301 Float_t IsolationVar; // isolation variable
302 Float_t IsolationVarRhoCorr; // isolation variable
303 Float_t SumPtCharged; // isolation variable
304 Float_t SumPtNeutral; // isolation variable
305 Float_t SumPtChargedPU; // isolation variable
306 Float_t SumPt; // isolation variable
307
308 Float_t D0; // track transverse impact parameter
309 Float_t DZ; // track longitudinal impact parameter
310 Float_t ErrorD0; // track transverse impact parameter error
311 Float_t ErrorDZ; // track longitudinal impact parameter error
312
313 static CompBase *fgCompare; //!
314 const CompBase *GetCompare() const { return fgCompare; }
315
316 TLorentzVector P4() const;
317
318 ClassDef(Electron, 4)
319};
320
321//---------------------------------------------------------------------------
322
323class Muon: public SortableObject
324{
325public:
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 Float_t D0; // track transverse impact parameter
344 Float_t DZ; // track longitudinal impact parameter
345 Float_t ErrorD0; // track transverse impact parameter error
346 Float_t ErrorDZ; // track longitudinal impact parameter error
347
348 static CompBase *fgCompare; //!
349 const CompBase *GetCompare() const { return fgCompare; }
350
351 TLorentzVector P4() const;
352
353 ClassDef(Muon, 4)
354};
355
356//---------------------------------------------------------------------------
357
358class Jet: public SortableObject
359{
360public:
361 Float_t PT; // jet transverse momentum
362 Float_t Eta; // jet pseudorapidity
363 Float_t Phi; // jet azimuthal angle
364
365 Float_t T; //particle arrival time of flight
366
367 Float_t Mass; // jet invariant mass
368
369 Float_t DeltaEta; // jet radius in pseudorapidity
370 Float_t DeltaPhi; // jet radius in azimuthal angle
371
372 UInt_t Flavor; // jet flavor
373 UInt_t FlavorAlgo; // jet flavor
374 UInt_t FlavorPhys; // jet flavor
375
376 UInt_t BTag; // 0 or 1 for a jet that has been tagged as containing a heavy quark
377 UInt_t BTagAlgo; // 0 or 1 for a jet that has been tagged as containing a heavy quark
378 UInt_t BTagPhys; // 0 or 1 for a jet that has been tagged as containing a heavy quark
379
380 UInt_t TauTag; // 0 or 1 for a jet that has been tagged as a tau
381 Float_t TauWeight; // probability for jet to be identified as tau
382
383 Int_t Charge; // tau charge
384
385 Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
386
387 Int_t NCharged; // number of charged constituents
388 Int_t NNeutrals; // number of neutral constituents
389
390 Float_t NeutralEnergyFraction; // charged energy fraction
391 Float_t ChargedEnergyFraction; // neutral energy fraction
392
393 Float_t Beta; // (sum pt of charged pile-up constituents)/(sum pt of charged constituents)
394 Float_t BetaStar; // (sum pt of charged constituents coming from hard interaction)/(sum pt of charged constituents)
395 Float_t MeanSqDeltaR; // average distance (squared) between constituent and jet weighted by pt (squared) of constituent
396 Float_t PTD; // average pt between constituent and jet weighted by pt of constituent
397 Float_t FracPt[5]; // (sum pt of constituents within a ring 0.1*i < DeltaR < 0.1*(i+1))/(sum pt of constituents)
398
399 Float_t Tau[5]; // N-subjettiness
400
401 TLorentzVector SoftDroppedJet;
402 TLorentzVector SoftDroppedSubJet1;
403 TLorentzVector SoftDroppedSubJet2;
404
405 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
406 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
407 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
408
409 Int_t NSubJetsTrimmed; // number of subjets trimmed
410 Int_t NSubJetsPruned; // number of subjets pruned
411 Int_t NSubJetsSoftDropped; // number of subjets soft-dropped
412
413 Double_t ExclYmerge23;
414 Double_t ExclYmerge34;
415 Double_t ExclYmerge45;
416 Double_t ExclYmerge56;
417
418 TRefArray Constituents; // references to constituents
419 TRefArray Particles; // references to generated particles
420
421 static CompBase *fgCompare; //!
422 const CompBase *GetCompare() const { return fgCompare; }
423
424 TLorentzVector P4() const;
425 TLorentzVector Area;
426
427 ClassDef(Jet, 4)
428};
429
430//---------------------------------------------------------------------------
431
432class Track: public SortableObject
433{
434public:
435 Int_t PID; // HEP ID number
436
437 Int_t Charge; // track charge
438
439 Float_t P; // track momentum
440 Float_t PT; // track transverse momentum
441 Float_t Eta; // track pseudorapidity
442 Float_t Phi; // track azimuthal angle
443 Float_t CtgTheta; // track cotangent of theta
444
445 Float_t EtaOuter; // track pseudorapidity at the tracker edge
446 Float_t PhiOuter; // track azimuthal angle at the tracker edge
447
448 Float_t T; // track vertex position (t component)
449 Float_t X; // track vertex position (x component)
450 Float_t Y; // track vertex position (y component)
451 Float_t Z; // track vertex position (z component)
452
453 Float_t TOuter; // track position (t component) at the tracker edge
454 Float_t XOuter; // track position (x component) at the tracker edge
455 Float_t YOuter; // track position (y component) at the tracker edge
456 Float_t ZOuter; // track position (z component) at the tracker edge
457
458 Float_t Td; // T coordinate of point of closest approach to vertex
459 Float_t Xd; // X coordinate of point of closest approach to vertex
460 Float_t Yd; // Y coordinate of point of closest approach to vertex
461 Float_t Zd; // Z coordinate of point of closest approach to vertex
462
463 Float_t TOFreco; // reconstructed time of flight
464 Float_t TOFgen; // gen time of flight
465
466 Float_t L; // track path length
467 Float_t D0; // track transverse impact parameter
468 Float_t DZ; // track longitudinal impact parameter
469
470 Float_t ErrorP; // track momentum error
471 Float_t ErrorPT; // track transverse momentum error
472 Float_t ErrorPhi; // track azimuthal angle error
473 Float_t ErrorCtgTheta; // track cotangent of theta error
474
475 Float_t ErrorT; // time measurement error
476 Float_t ErrorD0; // track transverse impact parameter error
477 Float_t ErrorDZ; // track longitudinal impact parameter error
478
479 TRef Particle; // reference to generated particle
480
481 Int_t VertexIndex; // reference to vertex
482
483 static CompBase *fgCompare; //!
484 const CompBase *GetCompare() const { return fgCompare; }
485
486 TLorentzVector P4() const;
487
488 ClassDef(Track, 4)
489};
490
491//---------------------------------------------------------------------------
492
493class Tower: public SortableObject
494{
495public:
496 Float_t ET; // calorimeter tower transverse energy
497 Float_t Eta; // calorimeter tower pseudorapidity
498 Float_t Phi; // calorimeter tower azimuthal angle
499
500 Float_t E; // calorimeter tower energy
501
502 Float_t T; // ecal deposit time, averaged by sqrt(EM energy) over all particles, not smeared
503 Int_t NTimeHits; // number of hits contributing to time measurement
504
505 Float_t Eem; // calorimeter tower electromagnetic energy
506 Float_t Ehad; // calorimeter tower hadronic energy
507
508 Float_t Edges[4]; // calorimeter tower edges
509
510 TRefArray Particles; // references to generated particles
511
512 static CompBase *fgCompare; //!
513 const CompBase *GetCompare() const { return fgCompare; }
514
515 TLorentzVector P4() const;
516
517 ClassDef(Tower, 2)
518};
519
520//---------------------------------------------------------------------------
521
522class ParticleFlowCandidate: public SortableObject
523{
524
525public:
526 Int_t PID; // HEP ID number
527
528 Int_t Charge; // track charge
529
530 Float_t E; // reconstructed energy [GeV]
531 Float_t P; // track momentum
532 Float_t PT; // track transverse momentum
533 Float_t Eta; // track pseudorapidity
534 Float_t Phi; // track azimuthal angle
535 Float_t CtgTheta; // track cotangent of theta
536
537 Float_t EtaOuter; // track pseudorapidity at the tracker edge
538 Float_t PhiOuter; // track azimuthal angle at the tracker edge
539
540 Float_t T; // track vertex position (t component)
541 Float_t X; // track vertex position (x component)
542 Float_t Y; // track vertex position (y component)
543 Float_t Z; // track vertex position (z component)
544
545 Float_t TOuter; // track position (t component) at the tracker edge
546 Float_t XOuter; // track position (x component) at the tracker edge
547 Float_t YOuter; // track position (y component) at the tracker edge
548 Float_t ZOuter; // track position (z component) at the tracker edge
549
550 Float_t Td; // X coordinate of point of closest approach to vertex
551 Float_t Xd; // X coordinate of point of closest approach to vertex
552 Float_t Yd; // Y coordinate of point of closest approach to vertex
553 Float_t Zd; // Z coordinate of point of closest approach to vertex
554
555 Float_t L; // track path length
556 Float_t D0; // track transverse impact parameter
557 Float_t DZ; // track longitudinal impact parameter
558
559 Float_t ErrorP; // track momentum error
560 Float_t ErrorPT; // track transverse momentum error
561 Float_t ErrorPhi; // track azimuthal angle error
562 Float_t ErrorCtgTheta; // track cotangent of theta error
563
564 Float_t ErrorT; // time measurement error
565 Float_t ErrorD0; // track transverse impact parameter error
566 Float_t ErrorDZ; // track longitudinal impact parameter error
567
568 Int_t VertexIndex; // reference to vertex
569
570 static CompBase *fgCompare; //!
571 const CompBase *GetCompare() const { return fgCompare; }
572
573 TLorentzVector P4() const;
574
575 Int_t NTimeHits; // number of hits contributing to time measurement
576
577 Float_t Eem; // calorimeter tower electromagnetic energy
578 Float_t Ehad; // calorimeter tower hadronic energy
579
580 Float_t Edges[4]; // calorimeter tower edges
581
582 TRefArray Particles; // references to generated particles
583
584 ClassDef(ParticleFlowCandidate, 1)
585
586};
587
588//---------------------------------------------------------------------------
589
590class HectorHit: public SortableObject
591{
592public:
593 Float_t E; // reconstructed energy [GeV]
594
595 Float_t Tx; // angle of the momentum in the horizontal (x,z) plane [urad]
596 Float_t Ty; // angle of the momentum in the verical (y,z) plane [urad]
597
598 Float_t T; // time of flight to the detector [s]
599
600 Float_t X; // horizontal distance to the beam [um]
601 Float_t Y; // vertical distance to the beam [um]
602 Float_t S; // distance to the interaction point [m]
603
604 TRef Particle; // reference to generated particle
605
606 static CompBase *fgCompare; //!
607 const CompBase *GetCompare() const { return fgCompare; }
608
609 ClassDef(HectorHit, 1)
610};
611
612//---------------------------------------------------------------------------
613
614class Candidate: public SortableObject
615{
616 friend class DelphesFactory;
617
618public:
619 Candidate();
620
621 Int_t PID;
622
623 Int_t Status;
624 Int_t M1, M2, D1, D2;
625
626 Int_t Charge;
627
628 Float_t Mass;
629
630 Int_t IsPU;
631 Int_t IsRecoPU;
632
633 Int_t IsConstituent;
634 Int_t IsFromConversion;
635
636 UInt_t Flavor;
637 UInt_t FlavorAlgo;
638 UInt_t FlavorPhys;
639
640 UInt_t BTag;
641 UInt_t BTagAlgo;
642 UInt_t BTagPhys;
643
644 UInt_t TauTag;
645 Float_t TauWeight;
646
647 Float_t Eem;
648 Float_t Ehad;
649
650 Float_t Edges[4];
651 Float_t DeltaEta;
652 Float_t DeltaPhi;
653
654 TLorentzVector Momentum, Position, InitialPosition, PositionError, Area;
655
656 Float_t L; // path length
657 Float_t ErrorT; // path length
658 Float_t D0;
659 Float_t ErrorD0;
660 Float_t DZ;
661 Float_t ErrorDZ;
662 Float_t P;
663 Float_t ErrorP;
664 Float_t PT;
665 Float_t ErrorPT;
666 Float_t CtgTheta;
667 Float_t ErrorCtgTheta;
668 Float_t Phi;
669 Float_t ErrorPhi;
670
671 Float_t Td;
672 Float_t Xd;
673 Float_t Yd;
674 Float_t Zd;
675
676 // Vertexing variables
677 Float_t VertexingWeight;
678
679 // tracking resolution
680 Float_t TrackResolution;
681
682 // PileUpJetID variables
683
684 Int_t NCharged;
685 Int_t NNeutrals;
686 Float_t Beta;
687 Float_t BetaStar;
688 Float_t MeanSqDeltaR;
689 Float_t PTD;
690 Float_t FracPt[5];
691 Float_t NeutralEnergyFraction; // charged energy fraction
692 Float_t ChargedEnergyFraction; // neutral energy fraction
693
694
695 // Timing information
696
697 Int_t NTimeHits;
698 std::vector<std::pair<Float_t, Float_t> > ECalEnergyTimePairs;
699
700 // Isolation variables
701
702 Float_t IsolationVar;
703 Float_t IsolationVarRhoCorr;
704 Float_t SumPtCharged;
705 Float_t SumPtNeutral;
706 Float_t SumPtChargedPU;
707 Float_t SumPt;
708
709 // vertex variables
710
711 Int_t ClusterIndex;
712 Int_t ClusterNDF;
713 Double_t ClusterSigma;
714 Double_t SumPT2;
715 Double_t BTVSumPT2;
716 Double_t GenDeltaZ;
717 Double_t GenSumPT2;
718
719 // N-subjettiness variables
720
721 Float_t Tau[5];
722
723 // Other Substructure variables
724
725 TLorentzVector SoftDroppedJet;
726 TLorentzVector SoftDroppedSubJet1;
727 TLorentzVector SoftDroppedSubJet2;
728
729 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
730 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
731 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
732
733 Int_t NSubJetsTrimmed; // number of subjets trimmed
734 Int_t NSubJetsPruned; // number of subjets pruned
735 Int_t NSubJetsSoftDropped; // number of subjets soft-dropped
736
737 // Exclusive clustering variables
738 Double_t ExclYmerge23;
739 Double_t ExclYmerge34;
740 Double_t ExclYmerge45;
741 Double_t ExclYmerge56;
742
743 static CompBase *fgCompare; //!
744 const CompBase *GetCompare() const { return fgCompare; }
745
746 void AddCandidate(Candidate *object);
747 TObjArray *GetCandidates();
748
749 Bool_t Overlaps(const Candidate *object) const;
750
751 virtual void Copy(TObject &object) const;
752 virtual TObject *Clone(const char *newname = "") const;
753 virtual void Clear(Option_t *option = "");
754
755private:
756 DelphesFactory *fFactory; //!
757 TObjArray *fArray; //!
758
759 void SetFactory(DelphesFactory *factory) { fFactory = factory; }
760
761 ClassDef(Candidate, 6)
762};
763
764#endif // DelphesClasses_h
Note: See TracBrowser for help on using the repository browser.