Fork me on GitHub

source: git/classes/DelphesClasses.h@ 3db5282

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

included timing information in Calorimeter

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