Fork me on GitHub

source: git/classes/DelphesClasses.h@ c61b5ce

Last change on this file since c61b5ce was c61b5ce, checked in by christinaw97 <christina.wang@…>, 3 years ago

cleaning

  • Property mode set to 100644
File size: 24.6 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 Float_t decayX;
161 Float_t decayY;
162 Float_t decayZ;
163 Float_t decayT;
164
165 static CompBase *fgCompare; //!
166 const CompBase *GetCompare() const { return fgCompare; }
167
168 TLorentzVector P4() const;
169
170 ClassDef(GenParticle, 2)
171};
172
173//---------------------------------------------------------------------------
174
175class Vertex: public SortableObject
176{
177public:
178 Float_t T; // vertex position (t component)
179 Float_t X; // vertex position (x component)
180 Float_t Y; // vertex position (y component)
181 Float_t Z; // vertex position (z component)
182
183 Double_t ErrorT; // vertex position error (t component)
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)
187
188 Int_t Index; // vertex index
189 Int_t NDF; // number of degrees of freedom
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
197
198 TRefArray Constituents; // references to constituents
199
200 static CompBase *fgCompare; //!
201 const CompBase *GetCompare() const { return fgCompare; }
202
203 ClassDef(Vertex, 3)
204};
205
206//---------------------------------------------------------------------------
207
208class MissingET: public TObject
209{
210public:
211 Float_t MET; // mising transverse energy
212 Float_t Eta; // mising energy pseudorapidity
213 Float_t Phi; // mising energy azimuthal angle
214
215 TLorentzVector P4() const;
216
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
232class Rho: public TObject
233{
234public:
235 Float_t Rho; // rho energy density
236 Float_t Edges[2]; // pseudorapidity range edges
237
238 ClassDef(Rho, 1)
239};
240
241//---------------------------------------------------------------------------
242
243class Weight: public TObject
244{
245public:
246 Float_t Weight; // weight for the event
247
248 ClassDef(Weight, 1)
249};
250
251//---------------------------------------------------------------------------
252
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
261
262 Float_t T; // particle arrival time of flight
263
264 Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
265
266 TRefArray Particles; // references to generated particles
267
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
274
275 Int_t Status; // 1: prompt -- 2: non prompt -- 3: fake
276
277 static CompBase *fgCompare; //!
278 const CompBase *GetCompare() const { return fgCompare; }
279
280 TLorentzVector P4() const;
281
282 ClassDef(Photon, 4)
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
293
294 Float_t T; // particle arrival time of flight
295
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
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
308
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
314 static CompBase *fgCompare; //!
315 const CompBase *GetCompare() const { return fgCompare; }
316
317 TLorentzVector P4() const;
318
319 ClassDef(Electron, 4)
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
331 Float_t T; // particle arrival time of flight
332
333 Int_t Charge; // muon charge
334
335 TRef Particle; // reference to generated particle
336
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
343
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
349 static CompBase *fgCompare; //!
350 const CompBase *GetCompare() const { return fgCompare; }
351
352 TLorentzVector P4() const;
353
354 ClassDef(Muon, 4)
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
366 Float_t T; //particle arrival time of flight
367
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
373 UInt_t Flavor; // jet flavor
374 UInt_t FlavorAlgo; // jet flavor
375 UInt_t FlavorPhys; // jet flavor
376
377 UInt_t BTag; // 0 or 1 for a jet that has been tagged as containing a heavy quark
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
380
381 UInt_t TauTag; // 0 or 1 for a jet that has been tagged as a tau
382 Float_t TauWeight; // probability for jet to be identified as tau
383
384 Int_t Charge; // tau charge
385
386 Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
387
388 Int_t NCharged; // number of charged constituents
389 Int_t NNeutrals; // number of neutral constituents
390
391 Float_t NeutralEnergyFraction; // charged energy fraction
392 Float_t ChargedEnergyFraction; // neutral energy fraction
393
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)
399
400 Float_t Tau[5]; // N-subjettiness
401
402 TLorentzVector SoftDroppedJet;
403 TLorentzVector SoftDroppedSubJet1;
404 TLorentzVector SoftDroppedSubJet2;
405
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
414 Double_t ExclYmerge23;
415 Double_t ExclYmerge34;
416 Double_t ExclYmerge45;
417 Double_t ExclYmerge56;
418
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
425 TLorentzVector P4() const;
426 TLorentzVector Area;
427
428 ClassDef(Jet, 4)
429};
430
431//---------------------------------------------------------------------------
432
433class Track: public SortableObject
434{
435public:
436 Int_t PID; // HEP ID number
437
438 Int_t Charge; // track charge
439
440 Float_t P; // track momentum
441 Float_t PT; // track transverse momentum
442 Float_t Eta; // track pseudorapidity
443 Float_t Phi; // track azimuthal angle
444 Float_t CtgTheta; // track cotangent of theta
445 Float_t C; // track curvature inverse
446 Float_t Mass; // particle mass
447
448 Float_t EtaOuter; // track pseudorapidity at the tracker edge
449 Float_t PhiOuter; // track azimuthal angle at the tracker edge
450
451 Float_t T; // track vertex position (t component)
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
456 Float_t TOuter; // track position (t component) at the tracker edge
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
460
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
464
465 Float_t L; // track path length
466 Float_t D0; // track transverse impact parameter
467 Float_t DZ; // track longitudinal impact parameter
468 Float_t Nclusters; // Number of ionization clusters
469 Float_t dNdx; // Number of ionization clusters
470
471 Float_t ErrorP; // track momentum error
472 Float_t ErrorPT; // track transverse momentum error
473 Float_t ErrorPhi; // track azimuthal angle error
474 Float_t ErrorCtgTheta; // track cotangent of theta error
475
476 Float_t ErrorT; // time measurement error
477 Float_t ErrorD0; // track transverse impact parameter error
478 Float_t ErrorDZ; // track longitudinal impact parameter error
479 Float_t ErrorC; // track curvature error
480
481 // track covariance off-diagonal terms
482 Float_t ErrorD0Phi;
483 Float_t ErrorD0C;
484 Float_t ErrorD0DZ;
485 Float_t ErrorD0CtgTheta;
486 Float_t ErrorPhiC;
487 Float_t ErrorPhiDZ;
488 Float_t ErrorPhiCtgTheta ;
489 Float_t ErrorCDZ;
490 Float_t ErrorCCtgTheta;
491 Float_t ErrorDZCtgTheta;
492
493 TRef Particle; // reference to generated particle
494
495 Int_t VertexIndex; // reference to vertex
496
497 static CompBase *fgCompare; //!
498 const CompBase *GetCompare() const { return fgCompare; }
499
500 TLorentzVector P4() const;
501 TMatrixDSym CovarianceMatrix() const;
502
503 ClassDef(Track, 3)
504};
505
506//---------------------------------------------------------------------------
507
508class Tower: public SortableObject
509{
510public:
511 Float_t ET; // calorimeter tower transverse energy
512 Float_t Eta; // calorimeter tower pseudorapidity
513 Float_t Phi; // calorimeter tower azimuthal angle
514
515 Float_t E; // calorimeter tower energy
516
517 Float_t T; // ecal deposit time, averaged by sqrt(EM energy) over all particles, not smeared
518 Int_t NTimeHits; // number of hits contributing to time measurement
519
520 Float_t Eem; // calorimeter tower electromagnetic energy
521 Float_t Ehad; // calorimeter tower hadronic energy
522 Float_t Etrk; // total charged energy hitting tower
523
524 Float_t Edges[4]; // calorimeter tower edges
525
526 TRefArray Particles; // references to generated particles
527
528 static CompBase *fgCompare; //!
529 const CompBase *GetCompare() const { return fgCompare; }
530
531 TLorentzVector P4() const;
532
533 ClassDef(Tower, 3)
534};
535
536//---------------------------------------------------------------------------
537
538class ParticleFlowCandidate: public SortableObject
539{
540
541public:
542 Int_t PID; // HEP ID number
543
544 Int_t Charge; // track charge
545
546 Float_t E; // reconstructed energy [GeV]
547 Float_t P; // track momentum
548 Float_t PT; // track transverse momentum
549 Float_t Eta; // track pseudorapidity
550 Float_t Phi; // track azimuthal angle
551 Float_t CtgTheta; // track cotangent of theta
552 Float_t C; // track curvature inverse
553 Float_t Mass; // particle mass
554
555 Float_t EtaOuter; // track pseudorapidity at the tracker edge
556 Float_t PhiOuter; // track azimuthal angle at the tracker edge
557
558 Float_t T; // track vertex position (t component)
559 Float_t X; // track vertex position (x component)
560 Float_t Y; // track vertex position (y component)
561 Float_t Z; // track vertex position (z component)
562
563 Float_t TOuter; // track position (t component) at the tracker edge
564 Float_t XOuter; // track position (x component) at the tracker edge
565 Float_t YOuter; // track position (y component) at the tracker edge
566 Float_t ZOuter; // track position (z component) at the tracker edge
567
568 Float_t Xd; // X coordinate of point of closest approach to vertex
569 Float_t Yd; // Y coordinate of point of closest approach to vertex
570 Float_t Zd; // Z coordinate of point of closest approach to vertex
571
572 Float_t L; // track path length
573 Float_t D0; // track transverse impact parameter
574 Float_t DZ; // track longitudinal impact parameter
575 Float_t Nclusters; // Number of ionization clusters
576 Float_t dNdx; // Number of ionization clusters
577
578 Float_t ErrorP; // track momentum error
579 Float_t ErrorPT; // track transverse momentum error
580 Float_t ErrorPhi; // track azimuthal angle error
581 Float_t ErrorCtgTheta; // track cotangent of theta error
582
583 Float_t ErrorT; // time measurement error
584 Float_t ErrorD0; // track transverse impact parameter error
585 Float_t ErrorDZ; // track longitudinal impact parameter error
586 Float_t ErrorC; // track curvature error
587
588 // track covariance off-diagonal terms
589 Float_t ErrorD0Phi;
590 Float_t ErrorD0C;
591 Float_t ErrorD0DZ;
592 Float_t ErrorD0CtgTheta;
593 Float_t ErrorPhiC;
594 Float_t ErrorPhiDZ;
595 Float_t ErrorPhiCtgTheta ;
596 Float_t ErrorCDZ;
597 Float_t ErrorCCtgTheta;
598 Float_t ErrorDZCtgTheta;
599
600 Int_t VertexIndex; // reference to vertex
601
602 static CompBase *fgCompare; //!
603 const CompBase *GetCompare() const { return fgCompare; }
604
605 TLorentzVector P4() const;
606 TMatrixDSym CovarianceMatrix() const;
607
608 Int_t NTimeHits; // number of hits contributing to time measurement
609
610 Float_t Eem; // calorimeter tower electromagnetic energy
611 Float_t Ehad; // calorimeter tower hadronic energy
612 Float_t Etrk; // total charged energy hitting tower
613
614 Float_t Edges[4]; // calorimeter tower edges
615
616 TRefArray Particles; // references to generated particles
617
618 ClassDef(ParticleFlowCandidate, 3)
619
620};
621
622//---------------------------------------------------------------------------
623
624class HectorHit: public SortableObject
625{
626public:
627 Float_t E; // reconstructed energy [GeV]
628
629 Float_t Tx; // angle of the momentum in the horizontal (x,z) plane [urad]
630 Float_t Ty; // angle of the momentum in the verical (y,z) plane [urad]
631
632 Float_t T; // time of flight to the detector [s]
633
634 Float_t X; // horizontal distance to the beam [um]
635 Float_t Y; // vertical distance to the beam [um]
636 Float_t S; // distance to the interaction point [m]
637
638 TRef Particle; // reference to generated particle
639
640 static CompBase *fgCompare; //!
641 const CompBase *GetCompare() const { return fgCompare; }
642
643 ClassDef(HectorHit, 1)
644};
645//---------------------------------------------------------------------------
646
647class CscCluster: public SortableObject
648{
649public:
650 Float_t Eta; // eta of LLP
651 Float_t Phi; // phi of LLP
652 Float_t PT; // pt of LLP
653 Float_t Px;// px of LLP
654 Float_t Py;// py of LLP
655 Float_t Pz;// pz of LLP
656 Float_t E; // E of LLP
657 Float_t Ehad; // had energy of LLP
658 Float_t Eem; // em energy of LLP
659 Float_t pid; // LLP pid
660 Float_t T; // LLP decay time-photon travel time
661 Float_t X; // LLP decay x
662 Float_t Y; // LLP decay y
663 Float_t Z; // LLP decay z
664 Float_t R; // LLP decay z
665 Float_t beta; // LLP beta
666 Float_t ctau; //LLP ctau
667
668
669 static CompBase *fgCompare; //!
670 const CompBase *GetCompare() const { return fgCompare; }
671
672 ClassDef(CscCluster, 4)
673};
674
675//---------------------------------------------------------------------------
676
677class Candidate: public SortableObject
678{
679 friend class DelphesFactory;
680
681public:
682 Candidate();
683
684 Int_t PID;
685
686 Int_t Status;
687 Int_t M1, M2, D1, D2;
688
689 Int_t Charge;
690
691 Float_t Mass;
692
693 Int_t IsPU;
694 Int_t IsRecoPU;
695
696 Int_t IsConstituent;
697 Int_t IsFromConversion;
698
699 UInt_t Flavor;
700 UInt_t FlavorAlgo;
701 UInt_t FlavorPhys;
702
703 UInt_t BTag;
704 UInt_t BTagAlgo;
705 UInt_t BTagPhys;
706
707 UInt_t TauTag;
708 Float_t TauWeight;
709
710 Float_t Eem;
711 Float_t Ehad;
712 Float_t Etrk;
713
714 Float_t Edges[4];
715 Float_t DeltaEta;
716 Float_t DeltaPhi;
717
718 TLorentzVector Momentum, Position, InitialPosition, DecayPosition, PositionError, Area;
719
720 Float_t L; // path length
721 Float_t DZ;
722 Float_t ErrorDZ;
723 Float_t ErrorT; // path length
724 Float_t D0;
725 Float_t ErrorD0;
726 Float_t C;
727 Float_t ErrorC;
728 Float_t P;
729 Float_t ErrorP;
730 Float_t PT;
731 Float_t ErrorPT;
732 Float_t CtgTheta;
733 Float_t ErrorCtgTheta;
734 Float_t Phi;
735 Float_t ErrorPhi;
736
737 Float_t Nclusters; // Number of ionization clusters
738 Float_t dNdx; // Number of ionization clusters per unit length
739
740 Float_t Xd;
741 Float_t Yd;
742 Float_t Zd;
743
744 // tracking resolution
745
746 Float_t TrackResolution;
747
748 // PileUpJetID variables
749
750 Int_t NCharged;
751 Int_t NNeutrals;
752 Float_t Beta;
753 Float_t BetaStar;
754 Float_t MeanSqDeltaR;
755 Float_t PTD;
756 Float_t FracPt[5];
757 Float_t NeutralEnergyFraction; // charged energy fraction
758 Float_t ChargedEnergyFraction; // neutral energy fraction
759
760
761 // Timing information
762
763 Int_t NTimeHits;
764 std::vector<std::pair<Float_t, Float_t> > ECalEnergyTimePairs;
765
766 // Isolation variables
767
768 Float_t IsolationVar;
769 Float_t IsolationVarRhoCorr;
770 Float_t SumPtCharged;
771 Float_t SumPtNeutral;
772 Float_t SumPtChargedPU;
773 Float_t SumPt;
774
775 // ACTS compliant 6x6 track covariance (D0, phi, Curvature, dz, ctg(theta))
776
777 TMatrixDSym TrackCovariance;
778
779 // vertex variables
780
781 Int_t ClusterIndex;
782 Int_t ClusterNDF;
783 Double_t ClusterSigma;
784 Double_t SumPT2;
785 Double_t BTVSumPT2;
786 Double_t GenDeltaZ;
787 Double_t GenSumPT2;
788
789 // N-subjettiness variables
790
791 Float_t Tau[5];
792
793 // Other Substructure variables
794
795 TLorentzVector SoftDroppedJet;
796 TLorentzVector SoftDroppedSubJet1;
797 TLorentzVector SoftDroppedSubJet2;
798
799 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
800 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
801 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
802
803 Int_t NSubJetsTrimmed; // number of subjets trimmed
804 Int_t NSubJetsPruned; // number of subjets pruned
805 Int_t NSubJetsSoftDropped; // number of subjets soft-dropped
806
807 // Exclusive clustering variables
808 Double_t ExclYmerge23;
809 Double_t ExclYmerge34;
810 Double_t ExclYmerge45;
811 Double_t ExclYmerge56;
812
813 // event characteristics variables
814 Double_t ParticleDensity; // particle multiplicity density in the proximity of the particle
815
816 static CompBase *fgCompare; //!
817 const CompBase *GetCompare() const { return fgCompare; }
818
819 void AddCandidate(Candidate *object);
820 TObjArray *GetCandidates();
821
822 Bool_t Overlaps(const Candidate *object) const;
823
824 virtual void Copy(TObject &object) const;
825 virtual TObject *Clone(const char *newname = "") const;
826 virtual void Clear(Option_t *option = "");
827
828private:
829 DelphesFactory *fFactory; //!
830 TObjArray *fArray; //!
831
832 void SetFactory(DelphesFactory *factory) { fFactory = factory; }
833
834 ClassDef(Candidate, 6)
835};
836
837#endif // DelphesClasses_h
Note: See TracBrowser for help on using the repository browser.