Fork me on GitHub

source: git/classes/DelphesClasses.h@ 29b722a

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

Merge branch 'master' of github.com:Christinaw97/delphes into HEAD

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