Fork me on GitHub

source: git/classes/DelphesClasses.h@ dd514ae

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

resolve conflict in DelphesClasses.h

  • Property mode set to 100644
File size: 16.1 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 "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;
52 Float_t ProcTime;
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
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
97class HepMCEvent: public Event
98{
99public:
100
101 Int_t ProcessID; // unique signal process id | signal_process_id()
102 Int_t MPI; // number of multi parton interactions | mpi ()
103
104 Float_t Weight; // weight for the event
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, 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
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 PT; // particle transverse momentum
150 Float_t Eta; // particle pseudorapidity
151 Float_t Phi; // particle azimuthal angle
152
153 Float_t Rapidity; // particle rapidity
154
155 Float_t T; // particle vertex position (t component) | hepevt.vhep[number][3]
156 Float_t X; // particle vertex position (x component) | hepevt.vhep[number][0]
157 Float_t Y; // particle vertex position (y component) | hepevt.vhep[number][1]
158 Float_t Z; // particle vertex position (z component) | hepevt.vhep[number][2]
159
160 static CompBase *fgCompare; //!
161 const CompBase *GetCompare() const { return fgCompare; }
162
163 TLorentzVector P4() const;
164
165 ClassDef(GenParticle, 1)
166};
167
168//---------------------------------------------------------------------------
169
170class Vertex: public TObject
171{
172public:
173 Float_t T; // vertex position (t component)
174 Float_t X; // vertex position (x component)
175 Float_t Y; // vertex position (y component)
176 Float_t Z; // vertex position (z component)
177
178 ClassDef(Vertex, 1)
179};
180
181//---------------------------------------------------------------------------
182
183class MissingET: public TObject
184{
185public:
186 Float_t MET; // mising transverse energy
187 Float_t Eta; // mising energy pseudorapidity
188 Float_t Phi; // mising energy azimuthal angle
189
190 TLorentzVector P4() const;
191
192 ClassDef(MissingET, 1)
193};
194
195//---------------------------------------------------------------------------
196
197class ScalarHT: public TObject
198{
199public:
200 Float_t HT; // scalar sum of transverse momenta
201
202 ClassDef(ScalarHT, 1)
203};
204
205//---------------------------------------------------------------------------
206
207class Rho: public TObject
208{
209public:
210 Float_t Rho; // rho energy density
211 Float_t Edges[2]; // pseudorapidity range edges
212
213 ClassDef(Rho, 1)
214};
215
216//---------------------------------------------------------------------------
217
218class Weight: public TObject
219{
220public:
221 Float_t Weight; // weight for the event
222
223 ClassDef(Weight, 1)
224};
225
226//---------------------------------------------------------------------------
227
228class Photon: public SortableObject
229{
230public:
231
232 Float_t PT; // photon transverse momentum
233 Float_t Eta; // photon pseudorapidity
234 Float_t Phi; // photon azimuthal angle
235
236 Float_t E; // photon energy
237
238 Float_t T; //particle arrival time of flight
239
240 Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
241
242 TRefArray Particles; // references to generated particles
243
244 // Isolation variables
245
246 Float_t IsolationVar;
247 Float_t IsolationVarRhoCorr;
248 Float_t SumPtCharged;
249 Float_t SumPtNeutral;
250 Float_t SumPtChargedPU;
251 Float_t SumPt;
252
253 static CompBase *fgCompare; //!
254 const CompBase *GetCompare() const { return fgCompare; }
255
256 TLorentzVector P4() const;
257
258 ClassDef(Photon, 2)
259};
260
261//---------------------------------------------------------------------------
262
263class Electron: public SortableObject
264{
265public:
266
267 Float_t PT; // electron transverse momentum
268 Float_t Eta; // electron pseudorapidity
269 Float_t Phi; // electron azimuthal angle
270
271 Float_t T; //particle arrival time of flight
272
273 Int_t Charge; // electron charge
274
275 Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
276
277 TRef Particle; // reference to generated particle
278
279 // Isolation variables
280
281 Float_t IsolationVar;
282 Float_t IsolationVarRhoCorr;
283 Float_t SumPtCharged;
284 Float_t SumPtNeutral;
285 Float_t SumPtChargedPU;
286 Float_t SumPt;
287
288 static CompBase *fgCompare; //!
289 const CompBase *GetCompare() const { return fgCompare; }
290
291 TLorentzVector P4() const;
292
293 ClassDef(Electron, 2)
294};
295
296//---------------------------------------------------------------------------
297
298class Muon: public SortableObject
299{
300public:
301
302 Float_t PT; // muon transverse momentum
303 Float_t Eta; // muon pseudorapidity
304 Float_t Phi; // muon azimuthal angle
305
306 Float_t T; //particle arrival time of flight
307
308 Int_t Charge; // muon charge
309
310 TRef Particle; // reference to generated particle
311
312 // Isolation variables
313
314 Float_t IsolationVar;
315 Float_t IsolationVarRhoCorr;
316 Float_t SumPtCharged;
317 Float_t SumPtNeutral;
318 Float_t SumPtChargedPU;
319 Float_t SumPt;
320
321 static CompBase *fgCompare; //!
322 const CompBase *GetCompare() const { return fgCompare; }
323
324 TLorentzVector P4() const;
325
326 ClassDef(Muon, 2)
327};
328
329//---------------------------------------------------------------------------
330
331class Jet: public SortableObject
332{
333public:
334
335 Float_t PT; // jet transverse momentum
336 Float_t Eta; // jet pseudorapidity
337 Float_t Phi; // jet azimuthal angle
338
339 Float_t T; //particle arrival time of flight
340
341 Float_t Mass; // jet invariant mass
342
343 Float_t DeltaEta; // jet radius in pseudorapidity
344 Float_t DeltaPhi; // jet radius in azimuthal angle
345
346 UInt_t BTag; // 0 or 1 for a jet that has been tagged as containing a heavy quark
347 UInt_t TauTag; // 0 or 1 for a jet that has been tagged as a tau
348
349 Int_t Charge; // tau charge
350
351 Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
352
353 Int_t NCharged; // number of charged constituents
354 Int_t NNeutrals; // number of neutral constituents
355 Float_t Beta; // (sum pt of charged pile-up constituents)/(sum pt of charged constituents)
356 Float_t BetaStar; // (sum pt of charged constituents coming from hard interaction)/(sum pt of charged constituents)
357 Float_t MeanSqDeltaR; // average distance (squared) between constituent and jet weighted by pt (squared) of constituent
358 Float_t PTD; // average pt between constituent and jet weighted by pt of constituent
359 Float_t FracPt[5]; // (sum pt of constituents within a ring 0.1*i < DeltaR < 0.1*(i+1))/(sum pt of constituents)
360
361 Float_t Tau[5]; // N-subjettiness
362
363 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
364 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
365 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
366
367 Int_t NSubJetsTrimmed; // number of subjets trimmed
368 Int_t NSubJetsPruned; // number of subjets pruned
369 Int_t NSubJetsSoftDropped; // number of subjets soft-dropped
370
371 TRefArray Constituents; // references to constituents
372 TRefArray Particles; // references to generated particles
373
374 static CompBase *fgCompare; //!
375 const CompBase *GetCompare() const { return fgCompare; }
376
377 TLorentzVector P4();
378 TLorentzVector Area;
379
380 ClassDef(Jet, 3)
381};
382
383//---------------------------------------------------------------------------
384
385class Track: public SortableObject
386{
387public:
388 Int_t PID; // HEP ID number
389
390 Int_t Charge; // track charge
391
392 Float_t PT; // track transverse momentum
393
394 Float_t Eta; // track pseudorapidity
395 Float_t Phi; // track azimuthal angle
396
397 Float_t EtaOuter; // track pseudorapidity at the tracker edge
398 Float_t PhiOuter; // track azimuthal angle at the tracker edge
399
400 Float_t X; // track vertex position (x component)
401 Float_t Y; // track vertex position (y component)
402 Float_t Z; // track vertex position (z component)
403 Float_t T; // track vertex position (z component)
404
405 Float_t XOuter; // track position (x component) at the tracker edge
406 Float_t YOuter; // track position (y component) at the tracker edge
407 Float_t ZOuter; // track position (z component) at the tracker edge
408 Float_t TOuter; // track position (z component) at the tracker edge
409
410 Float_t Dxy; // track signed transverse impact parameter
411 Float_t SDxy; // signed error on the track signed transverse impact parameter
412 Float_t Xd; // X coordinate of point of closest approach to vertex
413 Float_t Yd; // Y coordinate of point of closest approach to vertex
414 Float_t Zd; // Z coordinate of point of closest approach to vertex
415
416 TRef Particle; // reference to generated particle
417
418 static CompBase *fgCompare; //!
419 const CompBase *GetCompare() const { return fgCompare; }
420
421 TLorentzVector P4() const;
422
423 ClassDef(Track, 2)
424};
425
426//---------------------------------------------------------------------------
427
428class Tower: public SortableObject
429{
430public:
431 Float_t ET; // calorimeter tower transverse energy
432 Float_t Eta; // calorimeter tower pseudorapidity
433 Float_t Phi; // calorimeter tower azimuthal angle
434
435 Float_t E; // calorimeter tower energy
436
437 Float_t T; // ecal deposit time, averaged by sqrt(EM energy) over all particles, not smeared
438 Int_t Ntimes; // number of hits contributing to time measurement
439
440 Float_t Eem; // calorimeter tower electromagnetic energy
441 Float_t Ehad; // calorimeter tower hadronic energy
442
443 Float_t Edges[4]; // calorimeter tower edges
444
445 TRefArray Particles; // references to generated particles
446
447 static CompBase *fgCompare; //!
448 const CompBase *GetCompare() const { return fgCompare; }
449
450 TLorentzVector P4() const;
451
452 ClassDef(Tower, 1)
453};
454
455//---------------------------------------------------------------------------
456
457class HectorHit: public SortableObject
458{
459public:
460 Float_t E; // reconstructed energy [GeV]
461
462 Float_t Tx; // angle of the momentum in the horizontal (x,z) plane [urad]
463 Float_t Ty; // angle of the momentum in the verical (y,z) plane [urad]
464
465 Float_t T; // time of flight to the detector [s]
466
467 Float_t X; // horizontal distance to the beam [um]
468 Float_t Y; // vertical distance to the beam [um]
469 Float_t S; // distance to the interaction point [m]
470
471 TRef Particle; // reference to generated particle
472
473 static CompBase *fgCompare; //!
474 const CompBase *GetCompare() const { return fgCompare; }
475
476 ClassDef(HectorHit, 1)
477};
478
479//---------------------------------------------------------------------------
480
481class Candidate: public SortableObject
482{
483 friend class DelphesFactory;
484
485public:
486 Candidate();
487
488 Int_t PID;
489
490 Int_t Status;
491 Int_t M1, M2, D1, D2;
492
493 Int_t Charge;
494
495 Float_t Mass;
496
497 Int_t IsPU;
498 Int_t IsRecoPU;
499
500 Int_t IsConstituent;
501
502 UInt_t BTag;
503 UInt_t TauTag;
504
505 Float_t Eem;
506 Float_t Ehad;
507
508 Float_t Edges[4];
509 Float_t DeltaEta;
510 Float_t DeltaPhi;
511
512 TLorentzVector Momentum, Position, Area;
513
514 Float_t Dxy;
515 Float_t SDxy;
516 Float_t Xd;
517 Float_t Yd;
518 Float_t Zd;
519
520 // PileUpJetID variables
521
522 Int_t NCharged;
523 Int_t NNeutrals;
524 Float_t Beta;
525 Float_t BetaStar;
526 Float_t MeanSqDeltaR;
527 Float_t PTD;
528 Float_t FracPt[5];
529
530 //Timing information
531
532 Int_t Ntimes;
533 std::vector<std::pair<Float_t,Float_t> > Ecal_E_t;
534
535 // Isolation variables
536
537 Float_t IsolationVar;
538 Float_t IsolationVarRhoCorr;
539 Float_t SumPtCharged;
540 Float_t SumPtNeutral;
541 Float_t SumPtChargedPU;
542 Float_t SumPt;
543
544 // N-subjettiness variables
545
546 Float_t Tau[5];
547
548 // Other Substructure variables
549
550 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
551 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
552 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
553
554 Int_t NSubJetsTrimmed; // number of subjets trimmed
555 Int_t NSubJetsPruned; // number of subjets pruned
556 Int_t NSubJetsSoftDropped; // number of subjets soft-dropped
557
558
559 static CompBase *fgCompare; //!
560 const CompBase *GetCompare() const { return fgCompare; }
561
562 void AddCandidate(Candidate *object);
563 TObjArray *GetCandidates();
564
565 Bool_t Overlaps(const Candidate *object) const;
566
567 virtual void Copy(TObject &object) const;
568 virtual TObject *Clone(const char *newname = "") const;
569 virtual void Clear(Option_t* option = "");
570
571private:
572 DelphesFactory *fFactory; //!
573 TObjArray *fArray; //!
574
575 void SetFactory(DelphesFactory *factory) { fFactory = factory; }
576
577 ClassDef(Candidate, 3)
578};
579
580#endif // DelphesClasses_h
581
582
Note: See TracBrowser for help on using the repository browser.