Fork me on GitHub

source: git/classes/DelphesClasses.h@ f08474a

Last change on this file since f08474a was a95da74, checked in by michele <michele.selvaggi@…>, 4 years ago

added v0 of Ionisation cluster counting

  • Property mode set to 100644
File size: 23.5 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 "TMatrixDSym.h"
36#include "TObject.h"
37#include "TRef.h"
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
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 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
75 Float_t CrossSection; // cross-section (read from init, implemented only for Wizard evgen)
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
80 ClassDef(LHEFEvent, 3)
81};
82
83//---------------------------------------------------------------------------
84
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
96class HepMCEvent: public Event
97{
98public:
99 Int_t ProcessID; // unique signal process id | signal_process_id()
100 Int_t MPI; // number of multi parton interactions | mpi ()
101
102 Float_t Weight; // weight for the event
103 Float_t CrossSection; // cross-section in pb
104 Float_t CrossSectionError; // cross-section error in pb
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, 3)
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 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
160 static CompBase *fgCompare; //!
161 const CompBase *GetCompare() const { return fgCompare; }
162
163 TLorentzVector P4() const;
164
165 ClassDef(GenParticle, 2)
166};
167
168//---------------------------------------------------------------------------
169
170class Vertex: public SortableObject
171{
172public:
173 Float_t T; // vertex position (t component)
174 Float_t X; // vertex position (x component)
175 Float_t Y; // vertex position (y component)
176 Float_t Z; // vertex position (z component)
177
178 Double_t ErrorT; // vertex position error (t component)
179 Double_t ErrorX; // vertex position error (x component)
180 Double_t ErrorY; // vertex position error (y component)
181 Double_t ErrorZ; // vertex position error (z component)
182
183 Int_t Index; // vertex index
184 Int_t NDF; // number of degrees of freedom
185
186 Double_t Sigma; // vertex position (z component) error
187 Double_t SumPT2; // sum pt^2 of tracks attached to the vertex
188 Double_t GenSumPT2; // sum pt^2 of gen tracks attached to the vertex
189
190 Double_t GenDeltaZ; // distance in z to closest generated vertex
191 Double_t BTVSumPT2; // sum pt^2 of tracks attached to the secondary vertex
192
193 TRefArray Constituents; // references to constituents
194
195 static CompBase *fgCompare; //!
196 const CompBase *GetCompare() const { return fgCompare; }
197
198 ClassDef(Vertex, 3)
199};
200
201//---------------------------------------------------------------------------
202
203class MissingET: public TObject
204{
205public:
206 Float_t MET; // mising transverse energy
207 Float_t Eta; // mising energy pseudorapidity
208 Float_t Phi; // mising energy azimuthal angle
209
210 TLorentzVector P4() const;
211
212 ClassDef(MissingET, 1)
213};
214
215//---------------------------------------------------------------------------
216
217class ScalarHT: public TObject
218{
219public:
220 Float_t HT; // scalar sum of transverse momenta
221
222 ClassDef(ScalarHT, 1)
223};
224
225//---------------------------------------------------------------------------
226
227class Rho: public TObject
228{
229public:
230 Float_t Rho; // rho energy density
231 Float_t Edges[2]; // pseudorapidity range edges
232
233 ClassDef(Rho, 1)
234};
235
236//---------------------------------------------------------------------------
237
238class Weight: public TObject
239{
240public:
241 Float_t Weight; // weight for the event
242
243 ClassDef(Weight, 1)
244};
245
246//---------------------------------------------------------------------------
247
248class Photon: public SortableObject
249{
250public:
251 Float_t PT; // photon transverse momentum
252 Float_t Eta; // photon pseudorapidity
253 Float_t Phi; // photon azimuthal angle
254
255 Float_t E; // photon energy
256
257 Float_t T; // particle arrival time of flight
258
259 Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
260
261 TRefArray Particles; // references to generated particles
262
263 Float_t IsolationVar; // isolation variable
264 Float_t IsolationVarRhoCorr; // isolation variable
265 Float_t SumPtCharged; // isolation variable
266 Float_t SumPtNeutral; // isolation variable
267 Float_t SumPtChargedPU; // isolation variable
268 Float_t SumPt; // isolation variable
269
270 Int_t Status; // 1: prompt -- 2: non prompt -- 3: fake
271
272 static CompBase *fgCompare; //!
273 const CompBase *GetCompare() const { return fgCompare; }
274
275 TLorentzVector P4() const;
276
277 ClassDef(Photon, 4)
278};
279
280//---------------------------------------------------------------------------
281
282class Electron: public SortableObject
283{
284public:
285 Float_t PT; // electron transverse momentum
286 Float_t Eta; // electron pseudorapidity
287 Float_t Phi; // electron azimuthal angle
288
289 Float_t T; // particle arrival time of flight
290
291 Int_t Charge; // electron charge
292
293 Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
294
295 TRef Particle; // reference to generated particle
296
297 Float_t IsolationVar; // isolation variable
298 Float_t IsolationVarRhoCorr; // isolation variable
299 Float_t SumPtCharged; // isolation variable
300 Float_t SumPtNeutral; // isolation variable
301 Float_t SumPtChargedPU; // isolation variable
302 Float_t SumPt; // isolation variable
303
304 Float_t D0; // track transverse impact parameter
305 Float_t DZ; // track longitudinal impact parameter
306 Float_t ErrorD0; // track transverse impact parameter error
307 Float_t ErrorDZ; // track longitudinal impact parameter error
308
309 static CompBase *fgCompare; //!
310 const CompBase *GetCompare() const { return fgCompare; }
311
312 TLorentzVector P4() const;
313
314 ClassDef(Electron, 4)
315};
316
317//---------------------------------------------------------------------------
318
319class Muon: public SortableObject
320{
321public:
322 Float_t PT; // muon transverse momentum
323 Float_t Eta; // muon pseudorapidity
324 Float_t Phi; // muon azimuthal angle
325
326 Float_t T; // particle arrival time of flight
327
328 Int_t Charge; // muon charge
329
330 TRef Particle; // reference to generated particle
331
332 Float_t IsolationVar; // isolation variable
333 Float_t IsolationVarRhoCorr; // isolation variable
334 Float_t SumPtCharged; // isolation variable
335 Float_t SumPtNeutral; // isolation variable
336 Float_t SumPtChargedPU; // isolation variable
337 Float_t SumPt; // isolation variable
338
339 Float_t D0; // track transverse impact parameter
340 Float_t DZ; // track longitudinal impact parameter
341 Float_t ErrorD0; // track transverse impact parameter error
342 Float_t ErrorDZ; // track longitudinal impact parameter error
343
344 static CompBase *fgCompare; //!
345 const CompBase *GetCompare() const { return fgCompare; }
346
347 TLorentzVector P4() const;
348
349 ClassDef(Muon, 4)
350};
351
352//---------------------------------------------------------------------------
353
354class Jet: public SortableObject
355{
356public:
357 Float_t PT; // jet transverse momentum
358 Float_t Eta; // jet pseudorapidity
359 Float_t Phi; // jet azimuthal angle
360
361 Float_t T; //particle arrival time of flight
362
363 Float_t Mass; // jet invariant mass
364
365 Float_t DeltaEta; // jet radius in pseudorapidity
366 Float_t DeltaPhi; // jet radius in azimuthal angle
367
368 UInt_t Flavor; // jet flavor
369 UInt_t FlavorAlgo; // jet flavor
370 UInt_t FlavorPhys; // jet flavor
371
372 UInt_t BTag; // 0 or 1 for a jet that has been tagged as containing a heavy quark
373 UInt_t BTagAlgo; // 0 or 1 for a jet that has been tagged as containing a heavy quark
374 UInt_t BTagPhys; // 0 or 1 for a jet that has been tagged as containing a heavy quark
375
376 UInt_t TauTag; // 0 or 1 for a jet that has been tagged as a tau
377 Float_t TauWeight; // probability for jet to be identified as tau
378
379 Int_t Charge; // tau charge
380
381 Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
382
383 Int_t NCharged; // number of charged constituents
384 Int_t NNeutrals; // number of neutral constituents
385
386 Float_t NeutralEnergyFraction; // charged energy fraction
387 Float_t ChargedEnergyFraction; // neutral energy fraction
388
389 Float_t Beta; // (sum pt of charged pile-up constituents)/(sum pt of charged constituents)
390 Float_t BetaStar; // (sum pt of charged constituents coming from hard interaction)/(sum pt of charged constituents)
391 Float_t MeanSqDeltaR; // average distance (squared) between constituent and jet weighted by pt (squared) of constituent
392 Float_t PTD; // average pt between constituent and jet weighted by pt of constituent
393 Float_t FracPt[5]; // (sum pt of constituents within a ring 0.1*i < DeltaR < 0.1*(i+1))/(sum pt of constituents)
394
395 Float_t Tau[5]; // N-subjettiness
396
397 TLorentzVector SoftDroppedJet;
398 TLorentzVector SoftDroppedSubJet1;
399 TLorentzVector SoftDroppedSubJet2;
400
401 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
402 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
403 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
404
405 Int_t NSubJetsTrimmed; // number of subjets trimmed
406 Int_t NSubJetsPruned; // number of subjets pruned
407 Int_t NSubJetsSoftDropped; // number of subjets soft-dropped
408
409 Double_t ExclYmerge23;
410 Double_t ExclYmerge34;
411 Double_t ExclYmerge45;
412 Double_t ExclYmerge56;
413
414 TRefArray Constituents; // references to constituents
415 TRefArray Particles; // references to generated particles
416
417 static CompBase *fgCompare; //!
418 const CompBase *GetCompare() const { return fgCompare; }
419
420 TLorentzVector P4() const;
421 TLorentzVector Area;
422
423 ClassDef(Jet, 4)
424};
425
426//---------------------------------------------------------------------------
427
428class Track: public SortableObject
429{
430public:
431 Int_t PID; // HEP ID number
432
433 Int_t Charge; // track charge
434
435 Float_t P; // track momentum
436 Float_t PT; // track transverse momentum
437 Float_t Eta; // track pseudorapidity
438 Float_t Phi; // track azimuthal angle
439 Float_t CtgTheta; // track cotangent of theta
440 Float_t C; // track curvature inverse
441 Float_t Mass; // particle mass
442
443 Float_t EtaOuter; // track pseudorapidity at the tracker edge
444 Float_t PhiOuter; // track azimuthal angle at the tracker edge
445
446 Float_t T; // track vertex position (t component)
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
451 Float_t TOuter; // track position (t component) at the tracker edge
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
455
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
459
460 Float_t L; // track path length
461 Float_t D0; // track transverse impact parameter
462 Float_t DZ; // track longitudinal impact parameter
463 Float_t Nclusters; // Number of ionization clusters
464
465 Float_t ErrorP; // track momentum error
466 Float_t ErrorPT; // track transverse momentum error
467 Float_t ErrorPhi; // track azimuthal angle error
468 Float_t ErrorCtgTheta; // track cotangent of theta error
469
470 Float_t ErrorT; // time measurement error
471 Float_t ErrorD0; // track transverse impact parameter error
472 Float_t ErrorDZ; // track longitudinal impact parameter error
473 Float_t ErrorC; // track curvature error
474
475 // track covariance off-diagonal terms
476 Float_t ErrorD0Phi;
477 Float_t ErrorD0C;
478 Float_t ErrorD0DZ;
479 Float_t ErrorD0CtgTheta;
480 Float_t ErrorPhiC;
481 Float_t ErrorPhiDZ;
482 Float_t ErrorPhiCtgTheta ;
483 Float_t ErrorCDZ;
484 Float_t ErrorCCtgTheta;
485 Float_t ErrorDZCtgTheta;
486
487 TRef Particle; // reference to generated particle
488
489 Int_t VertexIndex; // reference to vertex
490
491 static CompBase *fgCompare; //!
492 const CompBase *GetCompare() const { return fgCompare; }
493
494 TLorentzVector P4() const;
495 TMatrixDSym CovarianceMatrix() const;
496
497 ClassDef(Track, 3)
498};
499
500//---------------------------------------------------------------------------
501
502class Tower: public SortableObject
503{
504public:
505 Float_t ET; // calorimeter tower transverse energy
506 Float_t Eta; // calorimeter tower pseudorapidity
507 Float_t Phi; // calorimeter tower azimuthal angle
508
509 Float_t E; // calorimeter tower energy
510
511 Float_t T; // ecal deposit time, averaged by sqrt(EM energy) over all particles, not smeared
512 Int_t NTimeHits; // number of hits contributing to time measurement
513
514 Float_t Eem; // calorimeter tower electromagnetic energy
515 Float_t Ehad; // calorimeter tower hadronic energy
516
517 Float_t Edges[4]; // calorimeter tower edges
518
519 TRefArray Particles; // references to generated particles
520
521 static CompBase *fgCompare; //!
522 const CompBase *GetCompare() const { return fgCompare; }
523
524 TLorentzVector P4() const;
525
526 ClassDef(Tower, 2)
527};
528
529//---------------------------------------------------------------------------
530
531class ParticleFlowCandidate: public SortableObject
532{
533
534public:
535 Int_t PID; // HEP ID number
536
537 Int_t Charge; // track charge
538
539 Float_t E; // reconstructed energy [GeV]
540 Float_t P; // track momentum
541 Float_t PT; // track transverse momentum
542 Float_t Eta; // track pseudorapidity
543 Float_t Phi; // track azimuthal angle
544 Float_t CtgTheta; // track cotangent of theta
545 Float_t C; // track curvature inverse
546 Float_t Mass; // particle mass
547
548 Float_t EtaOuter; // track pseudorapidity at the tracker edge
549 Float_t PhiOuter; // track azimuthal angle at the tracker edge
550
551 Float_t T; // track vertex position (t component)
552 Float_t X; // track vertex position (x component)
553 Float_t Y; // track vertex position (y component)
554 Float_t Z; // track vertex position (z component)
555
556 Float_t TOuter; // track position (t component) at the tracker edge
557 Float_t XOuter; // track position (x component) at the tracker edge
558 Float_t YOuter; // track position (y component) at the tracker edge
559 Float_t ZOuter; // track position (z component) at the tracker edge
560
561 Float_t Xd; // X coordinate of point of closest approach to vertex
562 Float_t Yd; // Y coordinate of point of closest approach to vertex
563 Float_t Zd; // Z coordinate of point of closest approach to vertex
564
565 Float_t L; // track path length
566 Float_t D0; // track transverse impact parameter
567 Float_t DZ; // track longitudinal impact parameter
568 Float_t Nclusters; // Number of ionization clusters
569
570 Float_t ErrorP; // track momentum error
571 Float_t ErrorPT; // track transverse momentum error
572 Float_t ErrorPhi; // track azimuthal angle error
573 Float_t ErrorCtgTheta; // track cotangent of theta error
574
575 Float_t ErrorT; // time measurement error
576 Float_t ErrorD0; // track transverse impact parameter error
577 Float_t ErrorDZ; // track longitudinal impact parameter error
578 Float_t ErrorC; // track curvature error
579
580 // track covariance off-diagonal terms
581 Float_t ErrorD0Phi;
582 Float_t ErrorD0C;
583 Float_t ErrorD0DZ;
584 Float_t ErrorD0CtgTheta;
585 Float_t ErrorPhiC;
586 Float_t ErrorPhiDZ;
587 Float_t ErrorPhiCtgTheta ;
588 Float_t ErrorCDZ;
589 Float_t ErrorCCtgTheta;
590 Float_t ErrorDZCtgTheta;
591
592 Int_t VertexIndex; // reference to vertex
593
594 static CompBase *fgCompare; //!
595 const CompBase *GetCompare() const { return fgCompare; }
596
597 TLorentzVector P4() const;
598 TMatrixDSym CovarianceMatrix() const;
599
600 Int_t NTimeHits; // number of hits contributing to time measurement
601
602 Float_t Eem; // calorimeter tower electromagnetic energy
603 Float_t Ehad; // calorimeter tower hadronic energy
604
605 Float_t Edges[4]; // calorimeter tower edges
606
607 TRefArray Particles; // references to generated particles
608
609 ClassDef(ParticleFlowCandidate, 2)
610
611};
612
613//---------------------------------------------------------------------------
614
615class HectorHit: public SortableObject
616{
617public:
618 Float_t E; // reconstructed energy [GeV]
619
620 Float_t Tx; // angle of the momentum in the horizontal (x,z) plane [urad]
621 Float_t Ty; // angle of the momentum in the verical (y,z) plane [urad]
622
623 Float_t T; // time of flight to the detector [s]
624
625 Float_t X; // horizontal distance to the beam [um]
626 Float_t Y; // vertical distance to the beam [um]
627 Float_t S; // distance to the interaction point [m]
628
629 TRef Particle; // reference to generated particle
630
631 static CompBase *fgCompare; //!
632 const CompBase *GetCompare() const { return fgCompare; }
633
634 ClassDef(HectorHit, 1)
635};
636
637//---------------------------------------------------------------------------
638
639class Candidate: public SortableObject
640{
641 friend class DelphesFactory;
642
643public:
644 Candidate();
645
646 Int_t PID;
647
648 Int_t Status;
649 Int_t M1, M2, D1, D2;
650
651 Int_t Charge;
652
653 Float_t Mass;
654
655 Int_t IsPU;
656 Int_t IsRecoPU;
657
658 Int_t IsConstituent;
659 Int_t IsFromConversion;
660
661 UInt_t Flavor;
662 UInt_t FlavorAlgo;
663 UInt_t FlavorPhys;
664
665 UInt_t BTag;
666 UInt_t BTagAlgo;
667 UInt_t BTagPhys;
668
669 UInt_t TauTag;
670 Float_t TauWeight;
671
672 Float_t Eem;
673 Float_t Ehad;
674
675 Float_t Edges[4];
676 Float_t DeltaEta;
677 Float_t DeltaPhi;
678
679 TLorentzVector Momentum, Position, InitialPosition, PositionError, Area;
680
681 Float_t L; // path length
682 Float_t DZ;
683 Float_t ErrorDZ;
684 Float_t ErrorT; // path length
685 Float_t D0;
686 Float_t ErrorD0;
687 Float_t C;
688 Float_t ErrorC;
689 Float_t P;
690 Float_t ErrorP;
691 Float_t PT;
692 Float_t ErrorPT;
693 Float_t CtgTheta;
694 Float_t ErrorCtgTheta;
695 Float_t Phi;
696 Float_t ErrorPhi;
697
698 Float_t Nclusters; // Number of ionization clusters
699
700 Float_t Xd;
701 Float_t Yd;
702 Float_t Zd;
703
704 // tracking resolution
705
706 Float_t TrackResolution;
707
708 // PileUpJetID variables
709
710 Int_t NCharged;
711 Int_t NNeutrals;
712 Float_t Beta;
713 Float_t BetaStar;
714 Float_t MeanSqDeltaR;
715 Float_t PTD;
716 Float_t FracPt[5];
717 Float_t NeutralEnergyFraction; // charged energy fraction
718 Float_t ChargedEnergyFraction; // neutral energy fraction
719
720
721 // Timing information
722
723 Int_t NTimeHits;
724 std::vector<std::pair<Float_t, Float_t> > ECalEnergyTimePairs;
725
726 // Isolation variables
727
728 Float_t IsolationVar;
729 Float_t IsolationVarRhoCorr;
730 Float_t SumPtCharged;
731 Float_t SumPtNeutral;
732 Float_t SumPtChargedPU;
733 Float_t SumPt;
734
735 // ACTS compliant 6x6 track covariance (D0, phi, Curvature, dz, ctg(theta))
736
737 TMatrixDSym TrackCovariance;
738
739 // vertex variables
740
741 Int_t ClusterIndex;
742 Int_t ClusterNDF;
743 Double_t ClusterSigma;
744 Double_t SumPT2;
745 Double_t BTVSumPT2;
746 Double_t GenDeltaZ;
747 Double_t GenSumPT2;
748
749 // N-subjettiness variables
750
751 Float_t Tau[5];
752
753 // Other Substructure variables
754
755 TLorentzVector SoftDroppedJet;
756 TLorentzVector SoftDroppedSubJet1;
757 TLorentzVector SoftDroppedSubJet2;
758
759 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
760 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
761 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
762
763 Int_t NSubJetsTrimmed; // number of subjets trimmed
764 Int_t NSubJetsPruned; // number of subjets pruned
765 Int_t NSubJetsSoftDropped; // number of subjets soft-dropped
766
767 // Exclusive clustering variables
768 Double_t ExclYmerge23;
769 Double_t ExclYmerge34;
770 Double_t ExclYmerge45;
771 Double_t ExclYmerge56;
772
773 // event characteristics variables
774 Double_t ParticleDensity; // particle multiplicity density in the proximity of the particle
775
776 static CompBase *fgCompare; //!
777 const CompBase *GetCompare() const { return fgCompare; }
778
779 void AddCandidate(Candidate *object);
780 TObjArray *GetCandidates();
781
782 Bool_t Overlaps(const Candidate *object) const;
783
784 virtual void Copy(TObject &object) const;
785 virtual TObject *Clone(const char *newname = "") const;
786 virtual void Clear(Option_t *option = "");
787
788private:
789 DelphesFactory *fFactory; //!
790 TObjArray *fArray; //!
791
792 void SetFactory(DelphesFactory *factory) { fFactory = factory; }
793
794 ClassDef(Candidate, 6)
795};
796
797#endif // DelphesClasses_h
Note: See TracBrowser for help on using the repository browser.