Fork me on GitHub

source: git/classes/DelphesClasses.h@ 5496767

ImprovedOutputFile Timing dual_readout llp
Last change on this file since 5496767 was 5496767, checked in by Michele Selvaggi <michele.selvaggi@…>, 8 years ago

added array of constituents to vertices

  • 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;
184
185 Int_t Index;
186 Int_t NDF;
187 Double_t Sigma;
188 Double_t SumPT2;
189 Double_t BTVSumPT2;
190 Double_t GenDeltaZ;
191 Double_t GenSumPT2;
192
[5496767]193 TRefArray Constituents; // references to constituents
194
[3c46e17]195 static CompBase *fgCompare; //!
196 const CompBase *GetCompare() const { return fgCompare; }
197
198 ClassDef(Vertex, 3)
[d07e957]199};
200
201//---------------------------------------------------------------------------
202
[d7d2da3]203class MissingET: public TObject
204{
205public:
206 Float_t MET; // mising transverse energy
[4ad7b96]207 Float_t Eta; // mising energy pseudorapidity
[d7d2da3]208 Float_t Phi; // mising energy azimuthal angle
209
[2b3ef28]210 TLorentzVector P4() const;
[4ad7b96]211
[d7d2da3]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
[71648c2]227class Rho: public TObject
228{
229public:
230 Float_t Rho; // rho energy density
[3b465ca]231 Float_t Edges[2]; // pseudorapidity range edges
[71648c2]232
233 ClassDef(Rho, 1)
234};
235
236//---------------------------------------------------------------------------
237
[2e229c9]238class Weight: public TObject
239{
240public:
241 Float_t Weight; // weight for the event
242
243 ClassDef(Weight, 1)
244};
245
246//---------------------------------------------------------------------------
247
[d7d2da3]248class Photon: public SortableObject
249{
250public:
251
252 Float_t PT; // photon transverse momentum
253 Float_t Eta; // photon pseudorapidity
254 Float_t Phi; // photon azimuthal angle
255
256 Float_t E; // photon energy
[da00c35]257
258 Float_t T; //particle arrival time of flight
259
[d7d2da3]260 Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
261
262 TRefArray Particles; // references to generated particles
263
[b62c2da]264 // Isolation variables
[edf10ba]265
[b62c2da]266 Float_t IsolationVar;
267 Float_t IsolationVarRhoCorr;
268 Float_t SumPtCharged;
269 Float_t SumPtNeutral;
270 Float_t SumPtChargedPU;
271 Float_t SumPt;
272
[d7d2da3]273 static CompBase *fgCompare; //!
274 const CompBase *GetCompare() const { return fgCompare; }
275
[2b3ef28]276 TLorentzVector P4() const;
[d7d2da3]277
[eb52a5d]278 ClassDef(Photon, 3)
[d7d2da3]279};
280
281//---------------------------------------------------------------------------
282
283class Electron: public SortableObject
284{
285public:
286
287 Float_t PT; // electron transverse momentum
288 Float_t Eta; // electron pseudorapidity
289 Float_t Phi; // electron azimuthal angle
[da00c35]290
291 Float_t T; //particle arrival time of flight
292
[d7d2da3]293 Int_t Charge; // electron charge
294
295 Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
296
297 TRef Particle; // reference to generated particle
298
[b62c2da]299 // Isolation variables
[edf10ba]300
[b62c2da]301 Float_t IsolationVar;
302 Float_t IsolationVarRhoCorr;
303 Float_t SumPtCharged;
304 Float_t SumPtNeutral;
305 Float_t SumPtChargedPU;
306 Float_t SumPt;
307
[d7d2da3]308 static CompBase *fgCompare; //!
309 const CompBase *GetCompare() const { return fgCompare; }
310
[2b3ef28]311 TLorentzVector P4() const;
[d7d2da3]312
[eb52a5d]313 ClassDef(Electron, 3)
[d7d2da3]314};
315
316//---------------------------------------------------------------------------
317
318class Muon: public SortableObject
319{
320public:
321
322 Float_t PT; // muon transverse momentum
323 Float_t Eta; // muon pseudorapidity
324 Float_t Phi; // muon azimuthal angle
325
[da00c35]326 Float_t T; //particle arrival time of flight
327
[d7d2da3]328 Int_t Charge; // muon charge
329
330 TRef Particle; // reference to generated particle
331
[b62c2da]332 // Isolation variables
[edf10ba]333
[b62c2da]334 Float_t IsolationVar;
335 Float_t IsolationVarRhoCorr;
336 Float_t SumPtCharged;
337 Float_t SumPtNeutral;
338 Float_t SumPtChargedPU;
339 Float_t SumPt;
340
[d7d2da3]341 static CompBase *fgCompare; //!
342 const CompBase *GetCompare() const { return fgCompare; }
343
[2b3ef28]344 TLorentzVector P4() const;
[d7d2da3]345
[eb52a5d]346 ClassDef(Muon, 3)
[d7d2da3]347};
348
349//---------------------------------------------------------------------------
350
351class Jet: public SortableObject
352{
353public:
354
355 Float_t PT; // jet transverse momentum
356 Float_t Eta; // jet pseudorapidity
357 Float_t Phi; // jet azimuthal angle
358
[da00c35]359 Float_t T; //particle arrival time of flight
360
[d7d2da3]361 Float_t Mass; // jet invariant mass
362
363 Float_t DeltaEta; // jet radius in pseudorapidity
364 Float_t DeltaPhi; // jet radius in azimuthal angle
365
[fe0273c]366 UInt_t Flavor;
367 UInt_t FlavorAlgo;
368 UInt_t FlavorPhys;
[edf10ba]369
[fe0273c]370 UInt_t BTag; // 0 or 1 for a jet that has been tagged as containing a heavy quark
[edf10ba]371 UInt_t BTagAlgo;
[fe0273c]372 UInt_t BTagPhys;
[edf10ba]373
[264bf40]374 UInt_t TauTag; // 0 or 1 for a jet that has been tagged as a tau
[d7d2da3]375
376 Int_t Charge; // tau charge
377
378 Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
379
[edf10ba]380 Int_t NCharged; // number of charged constituents
381 Int_t NNeutrals; // number of neutral constituents
382 Float_t Beta; // (sum pt of charged pile-up constituents)/(sum pt of charged constituents)
383 Float_t BetaStar; // (sum pt of charged constituents coming from hard interaction)/(sum pt of charged constituents)
384 Float_t MeanSqDeltaR; // average distance (squared) between constituent and jet weighted by pt (squared) of constituent
385 Float_t PTD; // average pt between constituent and jet weighted by pt of constituent
386 Float_t FracPt[5]; // (sum pt of constituents within a ring 0.1*i < DeltaR < 0.1*(i+1))/(sum pt of constituents)
[63178fb]387
[666d795]388 Float_t Tau[5]; // N-subjettiness
[edf10ba]389
390 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
391 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
392 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
393
394 Int_t NSubJetsTrimmed; // number of subjets trimmed
395 Int_t NSubJetsPruned; // number of subjets pruned
396 Int_t NSubJetsSoftDropped; // number of subjets soft-dropped
397
[e4c3fef]398 TRefArray Constituents; // references to constituents
399 TRefArray Particles; // references to generated particles
400
401 static CompBase *fgCompare; //!
402 const CompBase *GetCompare() const { return fgCompare; }
403
[8707eeb]404 TLorentzVector P4() const;
[ba1f1ee]405 TLorentzVector Area;
[24d005f]406
[666d795]407 ClassDef(Jet, 3)
[d7d2da3]408};
409
410//---------------------------------------------------------------------------
411
[3b465ca]412class Track: public SortableObject
[d7d2da3]413{
[3b465ca]414public:
[d7d2da3]415 Int_t PID; // HEP ID number
416
417 Int_t Charge; // track charge
418
419 Float_t Eta; // track pseudorapidity
[5496767]420
[d7d2da3]421 Float_t EtaOuter; // track pseudorapidity at the tracker edge
422 Float_t PhiOuter; // track azimuthal angle at the tracker edge
423
424 Float_t X; // track vertex position (x component)
425 Float_t Y; // track vertex position (y component)
426 Float_t Z; // track vertex position (z component)
[22dc7fd]427 Float_t T; // track vertex position (z component)
[d7d2da3]428
429 Float_t XOuter; // track position (x component) at the tracker edge
430 Float_t YOuter; // track position (y component) at the tracker edge
431 Float_t ZOuter; // track position (z component) at the tracker edge
[22dc7fd]432 Float_t TOuter; // track position (z component) at the tracker edge
[e4c3fef]433
[acd0621]434 Float_t L; // track path length
[5496767]435 Float_t ErrorT; // error on the time measurement
436
[80306e6]437 Float_t D0; // track signed transverse impact parameter
438 Float_t ErrorD0; // signed error on the track signed transverse impact parameter
[5496767]439
[acd0621]440 Float_t DZ; // track transverse momentum
441 Float_t ErrorDZ; // track transverse momentum error
[5496767]442
[acd0621]443 Float_t P; // track transverse momentum
444 Float_t ErrorP; // track transverse momentum error
[5496767]445
[acd0621]446 Float_t PT; // track transverse momentum
447 Float_t ErrorPT; // track transverse momentum error
[5496767]448
[acd0621]449 Float_t CtgTheta; // track transverse momentum
450 Float_t ErrorCtgTheta; // track transverse momentum error
[5496767]451
[acd0621]452 Float_t Phi; // track azimuthal angle
453 Float_t ErrorPhi; // track azimuthal angle
[5496767]454
[a0431dc]455 Float_t Xd; // X coordinate of point of closest approach to vertex
456 Float_t Yd; // Y coordinate of point of closest approach to vertex
457 Float_t Zd; // Z coordinate of point of closest approach to vertex
[d7d2da3]458
459 TRef Particle; // reference to generated particle
460
[2600216]461 Int_t VertexIndex; // reference to vertex
[5496767]462
[d7d2da3]463 static CompBase *fgCompare; //!
464 const CompBase *GetCompare() const { return fgCompare; }
465
[2b3ef28]466 TLorentzVector P4() const;
[d7d2da3]467
[e4c3fef]468 ClassDef(Track, 2)
[d7d2da3]469};
470
471//---------------------------------------------------------------------------
472
[3b465ca]473class Tower: public SortableObject
[d7d2da3]474{
475public:
476 Float_t ET; // calorimeter tower transverse energy
477 Float_t Eta; // calorimeter tower pseudorapidity
478 Float_t Phi; // calorimeter tower azimuthal angle
479
480 Float_t E; // calorimeter tower energy
481
[3db5282]482 Float_t T; // ecal deposit time, averaged by sqrt(EM energy) over all particles, not smeared
[839deb7]483 Int_t NTimeHits; // number of hits contributing to time measurement
[edf10ba]484
[d7d2da3]485 Float_t Eem; // calorimeter tower electromagnetic energy
486 Float_t Ehad; // calorimeter tower hadronic energy
487
488 Float_t Edges[4]; // calorimeter tower edges
489
490 TRefArray Particles; // references to generated particles
491
492 static CompBase *fgCompare; //!
493 const CompBase *GetCompare() const { return fgCompare; }
494
[2b3ef28]495 TLorentzVector P4() const;
[d7d2da3]496
[eb52a5d]497 ClassDef(Tower, 2)
[d7d2da3]498};
499
500//---------------------------------------------------------------------------
501
[8f7db23]502class HectorHit: public SortableObject
503{
504public:
505 Float_t E; // reconstructed energy [GeV]
506
507 Float_t Tx; // angle of the momentum in the horizontal (x,z) plane [urad]
508 Float_t Ty; // angle of the momentum in the verical (y,z) plane [urad]
509
510 Float_t T; // time of flight to the detector [s]
511
512 Float_t X; // horizontal distance to the beam [um]
513 Float_t Y; // vertical distance to the beam [um]
514 Float_t S; // distance to the interaction point [m]
515
[64a4950]516 TRef Particle; // reference to generated particle
517
[8f7db23]518 static CompBase *fgCompare; //!
519 const CompBase *GetCompare() const { return fgCompare; }
520
521 ClassDef(HectorHit, 1)
522};
523
524//---------------------------------------------------------------------------
525
[3b465ca]526class Candidate: public SortableObject
[d7d2da3]527{
528 friend class DelphesFactory;
529
530public:
531 Candidate();
532
533 Int_t PID;
534
535 Int_t Status;
536 Int_t M1, M2, D1, D2;
537
538 Int_t Charge;
539
540 Float_t Mass;
[3b465ca]541
[d7d2da3]542 Int_t IsPU;
[b62c2da]543 Int_t IsRecoPU;
[edf10ba]544
[d7d2da3]545 Int_t IsConstituent;
[839deb7]546
[5d2481f]547 Int_t IsFromConversion;
[839deb7]548
[fe0273c]549 UInt_t Flavor;
550 UInt_t FlavorAlgo;
551 UInt_t FlavorPhys;
[edf10ba]552
[fe0273c]553 UInt_t BTag;
[edf10ba]554 UInt_t BTagAlgo;
[fe0273c]555 UInt_t BTagPhys;
[edf10ba]556
[264bf40]557 UInt_t TauTag;
[d7d2da3]558
559 Float_t Eem;
560 Float_t Ehad;
561
562 Float_t Edges[4];
563 Float_t DeltaEta;
564 Float_t DeltaPhi;
565
[0e2f49b]566 TLorentzVector Momentum, Position, InitialPosition, PositionError, Area;
[80306e6]567
568 Float_t L; // path length
[28c722a]569 Float_t ErrorT; // path length
[80306e6]570 Float_t D0;
571 Float_t ErrorD0;
572 Float_t DZ;
573 Float_t ErrorDZ;
574 Float_t P;
575 Float_t ErrorP;
576 Float_t PT;
577 Float_t ErrorPT;
578 Float_t CtgTheta;
579 Float_t ErrorCtgTheta;
580 Float_t Phi;
581 Float_t ErrorPhi;
[e4c3fef]582
[839deb7]583 Float_t Xd;
584 Float_t Yd;
585 Float_t Zd;
[d7d2da3]586
[a98c7ef]587 // tracking resolution
[5496767]588
[a98c7ef]589 Float_t TrackResolution;
590
[da00c35]591 // PileUpJetID variables
592
[839deb7]593 Int_t NCharged;
594 Int_t NNeutrals;
595 Float_t Beta;
596 Float_t BetaStar;
597 Float_t MeanSqDeltaR;
598 Float_t PTD;
599 Float_t FracPt[5];
[edf10ba]600
[839deb7]601 // Timing information
[edf10ba]602
[839deb7]603 Int_t NTimeHits;
604 std::vector< std::pair< Float_t, Float_t > > ECalEnergyTimePairs;
[e4c3fef]605
[b62c2da]606 // Isolation variables
[edf10ba]607
[b62c2da]608 Float_t IsolationVar;
609 Float_t IsolationVarRhoCorr;
610 Float_t SumPtCharged;
611 Float_t SumPtNeutral;
612 Float_t SumPtChargedPU;
613 Float_t SumPt;
614
[0e2f49b]615 // vertex variables
[5496767]616
[0e2f49b]617 Int_t ClusterIndex;
618 Int_t ClusterNDF;
619 Double_t ClusterSigma;
620 Double_t SumPT2;
621 Double_t BTVSumPT2;
622 Double_t GenDeltaZ;
623 Double_t GenSumPT2;
624
[63178fb]625 // N-subjettiness variables
[e4c3fef]626
627 Float_t Tau[5];
[edf10ba]628
[de6d698]629 // Other Substructure variables
[edf10ba]630
631 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
632 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
633 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
634
635 Int_t NSubJetsTrimmed; // number of subjets trimmed
636 Int_t NSubJetsPruned; // number of subjets pruned
637 Int_t NSubJetsSoftDropped; // number of subjets soft-dropped
[de6d698]638
[e4c3fef]639
[d7d2da3]640 static CompBase *fgCompare; //!
641 const CompBase *GetCompare() const { return fgCompare; }
642
643 void AddCandidate(Candidate *object);
644 TObjArray *GetCandidates();
645
646 Bool_t Overlaps(const Candidate *object) const;
647
648 virtual void Copy(TObject &object) const;
649 virtual TObject *Clone(const char *newname = "") const;
[3b465ca]650 virtual void Clear(Option_t* option = "");
[d7d2da3]651
652private:
653 DelphesFactory *fFactory; //!
654 TObjArray *fArray; //!
[3b465ca]655
[d7d2da3]656 void SetFactory(DelphesFactory *factory) { fFactory = factory; }
657
[3c46e17]658 ClassDef(Candidate, 5)
[d7d2da3]659};
660
661#endif // DelphesClasses_h
662
663
Note: See TracBrowser for help on using the repository browser.