Fork me on GitHub

source: git/classes/DelphesClasses.h@ 21eab4f

Last change on this file since 21eab4f was 332025f, checked in by Michele Selvaggi <michele.selvaggi@…>, 8 years ago

add time error in Vertex class

  • Property mode set to 100644
File size: 17.6 KB
RevLine 
[b443089]1/*
2 * Delphes: a framework for fast simulation of a generic collider experiment
3 * Copyright (C) 2012-2014 Universite catholique de Louvain (UCL), Belgium
[1fa50c2]4 *
[b443089]5 * This program is free software: you can redistribute it and/or modify
6 * it under the terms of the GNU General Public License as published by
7 * the Free Software Foundation, either version 3 of the License, or
8 * (at your option) any later version.
[1fa50c2]9 *
[b443089]10 * This program is distributed in the hope that it will be useful,
11 * but WITHOUT ANY WARRANTY; without even the implied warranty of
12 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13 * GNU General Public License for more details.
[1fa50c2]14 *
[b443089]15 * You should have received a copy of the GNU General Public License
16 * along with this program. If not, see <http://www.gnu.org/licenses/>.
17 */
18
[d7d2da3]19#ifndef DelphesClasses_h
20#define DelphesClasses_h
21
22/**
23 *
24 * Definition of classes to be stored in the root tree.
25 * Function CompareXYZ sorts objects by the variable XYZ that MUST be
26 * present in the data members of the root tree class of the branch.
27 *
28 * \author P. Demin - UCL, Louvain-la-Neuve
29 *
30 */
31
32// Dependencies (#includes)
33
34#include "TRef.h"
35#include "TObject.h"
36#include "TRefArray.h"
37#include "TLorentzVector.h"
38
39#include "classes/SortableObject.h"
40
41class DelphesFactory;
42
43//---------------------------------------------------------------------------
44
45class Event: public TObject
46{
47public:
48
49 Long64_t Number; // event number
50
51 Float_t ReadTime;
[3b465ca]52 Float_t ProcTime;
[d7d2da3]53
54 ClassDef(Event, 1)
55};
56
57//---------------------------------------------------------------------------
58
59class LHCOEvent: public Event
60{
61public:
62
63 Int_t Trigger; // trigger word
64
65 ClassDef(LHCOEvent, 1)
66};
67
68//---------------------------------------------------------------------------
69
70class LHEFEvent: public Event
71{
72public:
73
74 Int_t ProcessID; // subprocess code for the event | hepup.IDPRUP
75
76 Float_t Weight; // weight for the event | hepup.XWGTUP
77 Float_t ScalePDF; // scale in GeV used in the calculation of the PDFs in the event | hepup.SCALUP
78 Float_t AlphaQED; // value of the QED coupling used in the event | hepup.AQEDUP
79 Float_t AlphaQCD; // value of the QCD coupling used in the event | hepup.AQCDUP
80
81 ClassDef(LHEFEvent, 2)
82};
83
84//---------------------------------------------------------------------------
85
[986d9d5]86class LHEFWeight: public TObject
87{
88public:
89 Int_t ID; // weight ID
90 Float_t Weight; // weight value
91
92 ClassDef(LHEFWeight, 1)
93};
94
95//---------------------------------------------------------------------------
96
[d7d2da3]97class HepMCEvent: public Event
98{
99public:
100
101 Int_t ProcessID; // unique signal process id | signal_process_id()
[3b465ca]102 Int_t MPI; // number of multi parton interactions | mpi ()
[d7d2da3]103
[59abd43]104 Float_t Weight; // weight for the event
105
[d7d2da3]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()
[3b465ca]111 Int_t ID2; // flavour code of second parton | pdf_info()->id2()
[d7d2da3]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, 2)
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
[3b465ca]133
[d7d2da3]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
[acd0621]149 Float_t D0;
150 Float_t DZ;
151 Float_t P;
152 Float_t PT;
153 Float_t CtgTheta;
154 Float_t Phi;
[d7d2da3]155 Float_t Eta; // particle pseudorapidity
156 Float_t Rapidity; // particle rapidity
157
158 Float_t T; // particle vertex position (t component) | hepevt.vhep[number][3]
159 Float_t X; // particle vertex position (x component) | hepevt.vhep[number][0]
160 Float_t Y; // particle vertex position (y component) | hepevt.vhep[number][1]
161 Float_t Z; // particle vertex position (z component) | hepevt.vhep[number][2]
162
163 static CompBase *fgCompare; //!
164 const CompBase *GetCompare() const { return fgCompare; }
[3b465ca]165
[2b3ef28]166 TLorentzVector P4() const;
[d7d2da3]167
168 ClassDef(GenParticle, 1)
169};
170
171//---------------------------------------------------------------------------
172
[f59a7b6]173class Vertex: public SortableObject
[d07e957]174{
175public:
176 Float_t T; // vertex position (t component)
177 Float_t X; // vertex position (x component)
178 Float_t Y; // vertex position (y component)
179 Float_t Z; // vertex position (z component)
[5496767]180
[0e2f49b]181 Double_t ErrorX;
182 Double_t ErrorY;
183 Double_t ErrorZ;
[332025f]184 Double_t ErrorT;
[0e2f49b]185
186 Int_t Index;
187 Int_t NDF;
188 Double_t Sigma;
189 Double_t SumPT2;
190 Double_t BTVSumPT2;
191 Double_t GenDeltaZ;
192 Double_t GenSumPT2;
193
[5496767]194 TRefArray Constituents; // references to constituents
195
[3c46e17]196 static CompBase *fgCompare; //!
197 const CompBase *GetCompare() const { return fgCompare; }
198
199 ClassDef(Vertex, 3)
[d07e957]200};
201
202//---------------------------------------------------------------------------
203
[d7d2da3]204class MissingET: public TObject
205{
206public:
207 Float_t MET; // mising transverse energy
[4ad7b96]208 Float_t Eta; // mising energy pseudorapidity
[d7d2da3]209 Float_t Phi; // mising energy azimuthal angle
210
[2b3ef28]211 TLorentzVector P4() const;
[4ad7b96]212
[d7d2da3]213 ClassDef(MissingET, 1)
214};
215
216//---------------------------------------------------------------------------
217
218class ScalarHT: public TObject
219{
220public:
221 Float_t HT; // scalar sum of transverse momenta
222
223 ClassDef(ScalarHT, 1)
224};
225
226//---------------------------------------------------------------------------
227
[71648c2]228class Rho: public TObject
229{
230public:
231 Float_t Rho; // rho energy density
[3b465ca]232 Float_t Edges[2]; // pseudorapidity range edges
[71648c2]233
234 ClassDef(Rho, 1)
235};
236
237//---------------------------------------------------------------------------
238
[2e229c9]239class Weight: public TObject
240{
241public:
242 Float_t Weight; // weight for the event
243
244 ClassDef(Weight, 1)
245};
246
247//---------------------------------------------------------------------------
248
[d7d2da3]249class Photon: public SortableObject
250{
251public:
252
253 Float_t PT; // photon transverse momentum
254 Float_t Eta; // photon pseudorapidity
255 Float_t Phi; // photon azimuthal angle
256
257 Float_t E; // photon energy
[da00c35]258
259 Float_t T; //particle arrival time of flight
260
[d7d2da3]261 Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
262
263 TRefArray Particles; // references to generated particles
264
[b62c2da]265 // Isolation variables
[edf10ba]266
[b62c2da]267 Float_t IsolationVar;
268 Float_t IsolationVarRhoCorr;
269 Float_t SumPtCharged;
270 Float_t SumPtNeutral;
271 Float_t SumPtChargedPU;
272 Float_t SumPt;
273
[d7d2da3]274 static CompBase *fgCompare; //!
275 const CompBase *GetCompare() const { return fgCompare; }
276
[2b3ef28]277 TLorentzVector P4() const;
[d7d2da3]278
[eb52a5d]279 ClassDef(Photon, 3)
[d7d2da3]280};
281
282//---------------------------------------------------------------------------
283
284class Electron: public SortableObject
285{
286public:
287
288 Float_t PT; // electron transverse momentum
289 Float_t Eta; // electron pseudorapidity
290 Float_t Phi; // electron azimuthal angle
[da00c35]291
292 Float_t T; //particle arrival time of flight
293
[d7d2da3]294 Int_t Charge; // electron charge
295
296 Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
297
298 TRef Particle; // reference to generated particle
299
[b62c2da]300 // Isolation variables
[edf10ba]301
[b62c2da]302 Float_t IsolationVar;
303 Float_t IsolationVarRhoCorr;
304 Float_t SumPtCharged;
305 Float_t SumPtNeutral;
306 Float_t SumPtChargedPU;
307 Float_t SumPt;
308
[d7d2da3]309 static CompBase *fgCompare; //!
310 const CompBase *GetCompare() const { return fgCompare; }
311
[2b3ef28]312 TLorentzVector P4() const;
[d7d2da3]313
[eb52a5d]314 ClassDef(Electron, 3)
[d7d2da3]315};
316
317//---------------------------------------------------------------------------
318
319class Muon: public SortableObject
320{
321public:
322
323 Float_t PT; // muon transverse momentum
324 Float_t Eta; // muon pseudorapidity
325 Float_t Phi; // muon azimuthal angle
326
[da00c35]327 Float_t T; //particle arrival time of flight
328
[d7d2da3]329 Int_t Charge; // muon charge
330
331 TRef Particle; // reference to generated particle
332
[b62c2da]333 // Isolation variables
[edf10ba]334
[b62c2da]335 Float_t IsolationVar;
336 Float_t IsolationVarRhoCorr;
337 Float_t SumPtCharged;
338 Float_t SumPtNeutral;
339 Float_t SumPtChargedPU;
340 Float_t SumPt;
341
[d7d2da3]342 static CompBase *fgCompare; //!
343 const CompBase *GetCompare() const { return fgCompare; }
344
[2b3ef28]345 TLorentzVector P4() const;
[d7d2da3]346
[eb52a5d]347 ClassDef(Muon, 3)
[d7d2da3]348};
349
350//---------------------------------------------------------------------------
351
352class Jet: public SortableObject
353{
354public:
355
356 Float_t PT; // jet transverse momentum
357 Float_t Eta; // jet pseudorapidity
358 Float_t Phi; // jet azimuthal angle
359
[da00c35]360 Float_t T; //particle arrival time of flight
361
[d7d2da3]362 Float_t Mass; // jet invariant mass
363
364 Float_t DeltaEta; // jet radius in pseudorapidity
365 Float_t DeltaPhi; // jet radius in azimuthal angle
366
[fe0273c]367 UInt_t Flavor;
368 UInt_t FlavorAlgo;
369 UInt_t FlavorPhys;
[edf10ba]370
[fe0273c]371 UInt_t BTag; // 0 or 1 for a jet that has been tagged as containing a heavy quark
[edf10ba]372 UInt_t BTagAlgo;
[fe0273c]373 UInt_t BTagPhys;
[edf10ba]374
[264bf40]375 UInt_t TauTag; // 0 or 1 for a jet that has been tagged as a tau
[d7d2da3]376
377 Int_t Charge; // tau charge
378
379 Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
380
[edf10ba]381 Int_t NCharged; // number of charged constituents
382 Int_t NNeutrals; // number of neutral constituents
383 Float_t Beta; // (sum pt of charged pile-up constituents)/(sum pt of charged constituents)
384 Float_t BetaStar; // (sum pt of charged constituents coming from hard interaction)/(sum pt of charged constituents)
385 Float_t MeanSqDeltaR; // average distance (squared) between constituent and jet weighted by pt (squared) of constituent
386 Float_t PTD; // average pt between constituent and jet weighted by pt of constituent
387 Float_t FracPt[5]; // (sum pt of constituents within a ring 0.1*i < DeltaR < 0.1*(i+1))/(sum pt of constituents)
[63178fb]388
[666d795]389 Float_t Tau[5]; // N-subjettiness
[edf10ba]390
391 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
392 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
393 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
394
395 Int_t NSubJetsTrimmed; // number of subjets trimmed
396 Int_t NSubJetsPruned; // number of subjets pruned
397 Int_t NSubJetsSoftDropped; // number of subjets soft-dropped
398
[e4c3fef]399 TRefArray Constituents; // references to constituents
400 TRefArray Particles; // references to generated particles
401
402 static CompBase *fgCompare; //!
403 const CompBase *GetCompare() const { return fgCompare; }
404
[8707eeb]405 TLorentzVector P4() const;
[ba1f1ee]406 TLorentzVector Area;
[24d005f]407
[666d795]408 ClassDef(Jet, 3)
[d7d2da3]409};
410
411//---------------------------------------------------------------------------
412
[3b465ca]413class Track: public SortableObject
[d7d2da3]414{
[3b465ca]415public:
[d7d2da3]416 Int_t PID; // HEP ID number
417
418 Int_t Charge; // track charge
419
420 Float_t Eta; // track pseudorapidity
[5496767]421
[d7d2da3]422 Float_t EtaOuter; // track pseudorapidity at the tracker edge
423 Float_t PhiOuter; // track azimuthal angle at the tracker edge
424
425 Float_t X; // track vertex position (x component)
426 Float_t Y; // track vertex position (y component)
427 Float_t Z; // track vertex position (z component)
[22dc7fd]428 Float_t T; // track vertex position (z component)
[d7d2da3]429
430 Float_t XOuter; // track position (x component) at the tracker edge
431 Float_t YOuter; // track position (y component) at the tracker edge
432 Float_t ZOuter; // track position (z component) at the tracker edge
[22dc7fd]433 Float_t TOuter; // track position (z component) at the tracker edge
[e4c3fef]434
[acd0621]435 Float_t L; // track path length
[5496767]436 Float_t ErrorT; // error on the time measurement
437
[80306e6]438 Float_t D0; // track signed transverse impact parameter
439 Float_t ErrorD0; // signed error on the track signed transverse impact parameter
[5496767]440
[acd0621]441 Float_t DZ; // track transverse momentum
442 Float_t ErrorDZ; // track transverse momentum error
[5496767]443
[acd0621]444 Float_t P; // track transverse momentum
445 Float_t ErrorP; // track transverse momentum error
[5496767]446
[acd0621]447 Float_t PT; // track transverse momentum
448 Float_t ErrorPT; // track transverse momentum error
[5496767]449
[acd0621]450 Float_t CtgTheta; // track transverse momentum
451 Float_t ErrorCtgTheta; // track transverse momentum error
[5496767]452
[acd0621]453 Float_t Phi; // track azimuthal angle
454 Float_t ErrorPhi; // track azimuthal angle
[5496767]455
[a0431dc]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
[d7d2da3]459
460 TRef Particle; // reference to generated particle
461
[2600216]462 Int_t VertexIndex; // reference to vertex
[5496767]463
[d7d2da3]464 static CompBase *fgCompare; //!
465 const CompBase *GetCompare() const { return fgCompare; }
466
[2b3ef28]467 TLorentzVector P4() const;
[d7d2da3]468
[e4c3fef]469 ClassDef(Track, 2)
[d7d2da3]470};
471
472//---------------------------------------------------------------------------
473
[3b465ca]474class Tower: public SortableObject
[d7d2da3]475{
476public:
477 Float_t ET; // calorimeter tower transverse energy
478 Float_t Eta; // calorimeter tower pseudorapidity
479 Float_t Phi; // calorimeter tower azimuthal angle
480
481 Float_t E; // calorimeter tower energy
482
[3db5282]483 Float_t T; // ecal deposit time, averaged by sqrt(EM energy) over all particles, not smeared
[839deb7]484 Int_t NTimeHits; // number of hits contributing to time measurement
[edf10ba]485
[d7d2da3]486 Float_t Eem; // calorimeter tower electromagnetic energy
487 Float_t Ehad; // calorimeter tower hadronic energy
488
489 Float_t Edges[4]; // calorimeter tower edges
490
491 TRefArray Particles; // references to generated particles
492
493 static CompBase *fgCompare; //!
494 const CompBase *GetCompare() const { return fgCompare; }
495
[2b3ef28]496 TLorentzVector P4() const;
[d7d2da3]497
[eb52a5d]498 ClassDef(Tower, 2)
[d7d2da3]499};
500
501//---------------------------------------------------------------------------
502
[8f7db23]503class HectorHit: public SortableObject
504{
505public:
506 Float_t E; // reconstructed energy [GeV]
507
508 Float_t Tx; // angle of the momentum in the horizontal (x,z) plane [urad]
509 Float_t Ty; // angle of the momentum in the verical (y,z) plane [urad]
510
511 Float_t T; // time of flight to the detector [s]
512
513 Float_t X; // horizontal distance to the beam [um]
514 Float_t Y; // vertical distance to the beam [um]
515 Float_t S; // distance to the interaction point [m]
516
[64a4950]517 TRef Particle; // reference to generated particle
518
[8f7db23]519 static CompBase *fgCompare; //!
520 const CompBase *GetCompare() const { return fgCompare; }
521
522 ClassDef(HectorHit, 1)
523};
524
525//---------------------------------------------------------------------------
526
[3b465ca]527class Candidate: public SortableObject
[d7d2da3]528{
529 friend class DelphesFactory;
530
531public:
532 Candidate();
533
534 Int_t PID;
535
536 Int_t Status;
537 Int_t M1, M2, D1, D2;
538
539 Int_t Charge;
540
541 Float_t Mass;
[3b465ca]542
[d7d2da3]543 Int_t IsPU;
[b62c2da]544 Int_t IsRecoPU;
[edf10ba]545
[d7d2da3]546 Int_t IsConstituent;
[839deb7]547
[5d2481f]548 Int_t IsFromConversion;
[839deb7]549
[fe0273c]550 UInt_t Flavor;
551 UInt_t FlavorAlgo;
552 UInt_t FlavorPhys;
[edf10ba]553
[fe0273c]554 UInt_t BTag;
[edf10ba]555 UInt_t BTagAlgo;
[fe0273c]556 UInt_t BTagPhys;
[edf10ba]557
[264bf40]558 UInt_t TauTag;
[d7d2da3]559
560 Float_t Eem;
561 Float_t Ehad;
562
563 Float_t Edges[4];
564 Float_t DeltaEta;
565 Float_t DeltaPhi;
566
[0e2f49b]567 TLorentzVector Momentum, Position, InitialPosition, PositionError, Area;
[80306e6]568
569 Float_t L; // path length
[28c722a]570 Float_t ErrorT; // path length
[80306e6]571 Float_t D0;
572 Float_t ErrorD0;
573 Float_t DZ;
574 Float_t ErrorDZ;
575 Float_t P;
576 Float_t ErrorP;
577 Float_t PT;
578 Float_t ErrorPT;
579 Float_t CtgTheta;
580 Float_t ErrorCtgTheta;
581 Float_t Phi;
582 Float_t ErrorPhi;
[e4c3fef]583
[839deb7]584 Float_t Xd;
585 Float_t Yd;
586 Float_t Zd;
[d7d2da3]587
[a98c7ef]588 // tracking resolution
[5496767]589
[a98c7ef]590 Float_t TrackResolution;
591
[da00c35]592 // PileUpJetID variables
593
[839deb7]594 Int_t NCharged;
595 Int_t NNeutrals;
596 Float_t Beta;
597 Float_t BetaStar;
598 Float_t MeanSqDeltaR;
599 Float_t PTD;
600 Float_t FracPt[5];
[edf10ba]601
[839deb7]602 // Timing information
[edf10ba]603
[839deb7]604 Int_t NTimeHits;
605 std::vector< std::pair< Float_t, Float_t > > ECalEnergyTimePairs;
[e4c3fef]606
[b62c2da]607 // Isolation variables
[edf10ba]608
[b62c2da]609 Float_t IsolationVar;
610 Float_t IsolationVarRhoCorr;
611 Float_t SumPtCharged;
612 Float_t SumPtNeutral;
613 Float_t SumPtChargedPU;
614 Float_t SumPt;
615
[0e2f49b]616 // vertex variables
[5496767]617
[0e2f49b]618 Int_t ClusterIndex;
619 Int_t ClusterNDF;
620 Double_t ClusterSigma;
621 Double_t SumPT2;
622 Double_t BTVSumPT2;
623 Double_t GenDeltaZ;
624 Double_t GenSumPT2;
625
[63178fb]626 // N-subjettiness variables
[e4c3fef]627
628 Float_t Tau[5];
[edf10ba]629
[de6d698]630 // Other Substructure variables
[edf10ba]631
632 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
633 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
634 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
635
636 Int_t NSubJetsTrimmed; // number of subjets trimmed
637 Int_t NSubJetsPruned; // number of subjets pruned
638 Int_t NSubJetsSoftDropped; // number of subjets soft-dropped
[de6d698]639
[e4c3fef]640
[d7d2da3]641 static CompBase *fgCompare; //!
642 const CompBase *GetCompare() const { return fgCompare; }
643
644 void AddCandidate(Candidate *object);
645 TObjArray *GetCandidates();
646
647 Bool_t Overlaps(const Candidate *object) const;
648
649 virtual void Copy(TObject &object) const;
650 virtual TObject *Clone(const char *newname = "") const;
[3b465ca]651 virtual void Clear(Option_t* option = "");
[d7d2da3]652
653private:
654 DelphesFactory *fFactory; //!
655 TObjArray *fArray; //!
[3b465ca]656
[d7d2da3]657 void SetFactory(DelphesFactory *factory) { fFactory = factory; }
658
[3c46e17]659 ClassDef(Candidate, 5)
[d7d2da3]660};
661
662#endif // DelphesClasses_h
663
664
Note: See TracBrowser for help on using the repository browser.