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 |
|
---|
41 | class DelphesFactory;
|
---|
42 |
|
---|
43 | //---------------------------------------------------------------------------
|
---|
44 |
|
---|
45 | class Event: public TObject
|
---|
46 | {
|
---|
47 | public:
|
---|
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 |
|
---|
59 | class LHCOEvent: public Event
|
---|
60 | {
|
---|
61 | public:
|
---|
62 |
|
---|
63 | Int_t Trigger; // trigger word
|
---|
64 |
|
---|
65 | ClassDef(LHCOEvent, 1)
|
---|
66 | };
|
---|
67 |
|
---|
68 | //---------------------------------------------------------------------------
|
---|
69 |
|
---|
70 | class LHEFEvent: public Event
|
---|
71 | {
|
---|
72 | public:
|
---|
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 |
|
---|
86 | class HepMCEvent: public Event
|
---|
87 | {
|
---|
88 | public:
|
---|
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 |
|
---|
115 | class GenParticle: public SortableObject
|
---|
116 | {
|
---|
117 | public:
|
---|
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() const;
|
---|
153 |
|
---|
154 | ClassDef(GenParticle, 1)
|
---|
155 | };
|
---|
156 |
|
---|
157 | //---------------------------------------------------------------------------
|
---|
158 |
|
---|
159 | class Vertex: public TObject
|
---|
160 | {
|
---|
161 | public:
|
---|
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 |
|
---|
172 | class MissingET: public TObject
|
---|
173 | {
|
---|
174 | public:
|
---|
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() const;
|
---|
180 |
|
---|
181 | ClassDef(MissingET, 1)
|
---|
182 | };
|
---|
183 |
|
---|
184 | //---------------------------------------------------------------------------
|
---|
185 |
|
---|
186 | class ScalarHT: public TObject
|
---|
187 | {
|
---|
188 | public:
|
---|
189 | Float_t HT; // scalar sum of transverse momenta
|
---|
190 |
|
---|
191 | ClassDef(ScalarHT, 1)
|
---|
192 | };
|
---|
193 |
|
---|
194 | //---------------------------------------------------------------------------
|
---|
195 |
|
---|
196 | class Rho: public TObject
|
---|
197 | {
|
---|
198 | public:
|
---|
199 | Float_t Rho; // rho energy density
|
---|
200 | Float_t Edges[2]; // pseudorapidity range edges
|
---|
201 |
|
---|
202 | ClassDef(Rho, 1)
|
---|
203 | };
|
---|
204 |
|
---|
205 | //---------------------------------------------------------------------------
|
---|
206 |
|
---|
207 | class Weight: public TObject
|
---|
208 | {
|
---|
209 | public:
|
---|
210 | Float_t Weight; // weight for the event
|
---|
211 |
|
---|
212 | ClassDef(Weight, 1)
|
---|
213 | };
|
---|
214 |
|
---|
215 | //---------------------------------------------------------------------------
|
---|
216 |
|
---|
217 | class Photon: public SortableObject
|
---|
218 | {
|
---|
219 | public:
|
---|
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 | static CompBase *fgCompare; //!
|
---|
234 | const CompBase *GetCompare() const { return fgCompare; }
|
---|
235 |
|
---|
236 | TLorentzVector P4() const;
|
---|
237 |
|
---|
238 | ClassDef(Photon, 2)
|
---|
239 | };
|
---|
240 |
|
---|
241 | //---------------------------------------------------------------------------
|
---|
242 |
|
---|
243 | class Electron: public SortableObject
|
---|
244 | {
|
---|
245 | public:
|
---|
246 |
|
---|
247 | Float_t PT; // electron transverse momentum
|
---|
248 | Float_t Eta; // electron pseudorapidity
|
---|
249 | Float_t Phi; // electron azimuthal angle
|
---|
250 |
|
---|
251 | Float_t T; //particle arrival time of flight
|
---|
252 |
|
---|
253 | Int_t Charge; // electron charge
|
---|
254 |
|
---|
255 | Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
|
---|
256 |
|
---|
257 | TRef Particle; // reference to generated particle
|
---|
258 |
|
---|
259 | static CompBase *fgCompare; //!
|
---|
260 | const CompBase *GetCompare() const { return fgCompare; }
|
---|
261 |
|
---|
262 | TLorentzVector P4() const;
|
---|
263 |
|
---|
264 | ClassDef(Electron, 2)
|
---|
265 | };
|
---|
266 |
|
---|
267 | //---------------------------------------------------------------------------
|
---|
268 |
|
---|
269 | class Muon: public SortableObject
|
---|
270 | {
|
---|
271 | public:
|
---|
272 |
|
---|
273 | Float_t PT; // muon transverse momentum
|
---|
274 | Float_t Eta; // muon pseudorapidity
|
---|
275 | Float_t Phi; // muon azimuthal angle
|
---|
276 |
|
---|
277 | Float_t T; //particle arrival time of flight
|
---|
278 |
|
---|
279 | Int_t Charge; // muon charge
|
---|
280 |
|
---|
281 | TRef Particle; // reference to generated particle
|
---|
282 |
|
---|
283 | static CompBase *fgCompare; //!
|
---|
284 | const CompBase *GetCompare() const { return fgCompare; }
|
---|
285 |
|
---|
286 | TLorentzVector P4() const;
|
---|
287 |
|
---|
288 | ClassDef(Muon, 2)
|
---|
289 | };
|
---|
290 |
|
---|
291 | //---------------------------------------------------------------------------
|
---|
292 |
|
---|
293 | class Jet: public SortableObject
|
---|
294 | {
|
---|
295 | public:
|
---|
296 |
|
---|
297 | Float_t PT; // jet transverse momentum
|
---|
298 | Float_t Eta; // jet pseudorapidity
|
---|
299 | Float_t Phi; // jet azimuthal angle
|
---|
300 |
|
---|
301 | Float_t T; //particle arrival time of flight
|
---|
302 |
|
---|
303 | Float_t Mass; // jet invariant mass
|
---|
304 |
|
---|
305 | Float_t DeltaEta; // jet radius in pseudorapidity
|
---|
306 | Float_t DeltaPhi; // jet radius in azimuthal angle
|
---|
307 |
|
---|
308 | UInt_t BTag; // 0 or 1 for a jet that has been tagged as containing a heavy quark
|
---|
309 | UInt_t TauTag; // 0 or 1 for a jet that has been tagged as a tau
|
---|
310 |
|
---|
311 | Int_t Charge; // tau charge
|
---|
312 |
|
---|
313 | Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
|
---|
314 |
|
---|
315 | Int_t NCharged; // number of charged constituents
|
---|
316 | Int_t NNeutrals; // number of neutral constituents
|
---|
317 | Float_t Beta; // (sum pt of charged pile-up constituents)/(sum pt of charged constituents)
|
---|
318 | Float_t BetaStar; // (sum pt of charged constituents coming from hard interaction)/(sum pt of charged constituents)
|
---|
319 | Float_t MeanSqDeltaR; // average distance (squared) between constituent and jet weighted by pt (squared) of constituent
|
---|
320 | Float_t PTD; // average pt between constituent and jet weighted by pt of constituent
|
---|
321 | Float_t FracPt[5]; // (sum pt of constituents within a ring 0.1*i < DeltaR < 0.1*(i+1))/(sum pt of constituents)
|
---|
322 |
|
---|
323 | Float_t Tau1; // 1-subjettiness
|
---|
324 | Float_t Tau2; // 2-subjettiness
|
---|
325 | Float_t Tau3; // 3-subjettiness
|
---|
326 | Float_t Tau4; // 4-subjettiness
|
---|
327 | Float_t Tau5; // 5-subjettiness
|
---|
328 |
|
---|
329 | TRefArray Constituents; // references to constituents
|
---|
330 | TRefArray Particles; // references to generated particles
|
---|
331 |
|
---|
332 | static CompBase *fgCompare; //!
|
---|
333 | const CompBase *GetCompare() const { return fgCompare; }
|
---|
334 |
|
---|
335 | TLorentzVector P4() const;
|
---|
336 |
|
---|
337 | ClassDef(Jet, 2)
|
---|
338 | };
|
---|
339 |
|
---|
340 | //---------------------------------------------------------------------------
|
---|
341 |
|
---|
342 | class Track: public SortableObject
|
---|
343 | {
|
---|
344 | public:
|
---|
345 | Int_t PID; // HEP ID number
|
---|
346 |
|
---|
347 | Int_t Charge; // track charge
|
---|
348 |
|
---|
349 | Float_t PT; // track transverse momentum
|
---|
350 |
|
---|
351 | Float_t Eta; // track pseudorapidity
|
---|
352 | Float_t Phi; // track azimuthal angle
|
---|
353 |
|
---|
354 | Float_t EtaOuter; // track pseudorapidity at the tracker edge
|
---|
355 | Float_t PhiOuter; // track azimuthal angle at the tracker edge
|
---|
356 |
|
---|
357 | Float_t X; // track vertex position (x component)
|
---|
358 | Float_t Y; // track vertex position (y component)
|
---|
359 | Float_t Z; // track vertex position (z component)
|
---|
360 | Float_t T; // track vertex position (z component)
|
---|
361 |
|
---|
362 | Float_t XOuter; // track position (x component) at the tracker edge
|
---|
363 | Float_t YOuter; // track position (y component) at the tracker edge
|
---|
364 | Float_t ZOuter; // track position (z component) at the tracker edge
|
---|
365 | Float_t TOuter; // track position (z component) at the tracker edge
|
---|
366 |
|
---|
367 | Float_t Dxy; // track signed transverse impact parameter
|
---|
368 | Float_t SDxy; // signed error on the track signed transverse impact parameter
|
---|
369 | Float_t Xd; // X coordinate of point of closest approach to vertex
|
---|
370 | Float_t Yd; // Y coordinate of point of closest approach to vertex
|
---|
371 | Float_t Zd; // Z coordinate of point of closest approach to vertex
|
---|
372 |
|
---|
373 | TRef Particle; // reference to generated particle
|
---|
374 |
|
---|
375 | static CompBase *fgCompare; //!
|
---|
376 | const CompBase *GetCompare() const { return fgCompare; }
|
---|
377 |
|
---|
378 | TLorentzVector P4() const;
|
---|
379 |
|
---|
380 | ClassDef(Track, 2)
|
---|
381 | };
|
---|
382 |
|
---|
383 | //---------------------------------------------------------------------------
|
---|
384 |
|
---|
385 | class Tower: public SortableObject
|
---|
386 | {
|
---|
387 | public:
|
---|
388 | Float_t ET; // calorimeter tower transverse energy
|
---|
389 | Float_t Eta; // calorimeter tower pseudorapidity
|
---|
390 | Float_t Phi; // calorimeter tower azimuthal angle
|
---|
391 |
|
---|
392 | Float_t E; // calorimeter tower energy
|
---|
393 |
|
---|
394 | Float_t T; //particle arrival time of flight
|
---|
395 |
|
---|
396 | Float_t Eem; // calorimeter tower electromagnetic energy
|
---|
397 | Float_t Ehad; // calorimeter tower hadronic energy
|
---|
398 |
|
---|
399 | Float_t Edges[4]; // calorimeter tower edges
|
---|
400 |
|
---|
401 | TRefArray Particles; // references to generated particles
|
---|
402 |
|
---|
403 | static CompBase *fgCompare; //!
|
---|
404 | const CompBase *GetCompare() const { return fgCompare; }
|
---|
405 |
|
---|
406 | TLorentzVector P4() const;
|
---|
407 |
|
---|
408 | ClassDef(Tower, 1)
|
---|
409 | };
|
---|
410 |
|
---|
411 | //---------------------------------------------------------------------------
|
---|
412 |
|
---|
413 | class HectorHit: public SortableObject
|
---|
414 | {
|
---|
415 | public:
|
---|
416 | Float_t E; // reconstructed energy [GeV]
|
---|
417 |
|
---|
418 | Float_t Tx; // angle of the momentum in the horizontal (x,z) plane [urad]
|
---|
419 | Float_t Ty; // angle of the momentum in the verical (y,z) plane [urad]
|
---|
420 |
|
---|
421 | Float_t T; // time of flight to the detector [s]
|
---|
422 |
|
---|
423 | Float_t X; // horizontal distance to the beam [um]
|
---|
424 | Float_t Y; // vertical distance to the beam [um]
|
---|
425 | Float_t S; // distance to the interaction point [m]
|
---|
426 |
|
---|
427 | TRef Particle; // reference to generated particle
|
---|
428 |
|
---|
429 | static CompBase *fgCompare; //!
|
---|
430 | const CompBase *GetCompare() const { return fgCompare; }
|
---|
431 |
|
---|
432 | ClassDef(HectorHit, 1)
|
---|
433 | };
|
---|
434 |
|
---|
435 | //---------------------------------------------------------------------------
|
---|
436 |
|
---|
437 | class Candidate: public SortableObject
|
---|
438 | {
|
---|
439 | friend class DelphesFactory;
|
---|
440 |
|
---|
441 | public:
|
---|
442 | Candidate();
|
---|
443 |
|
---|
444 | Int_t PID;
|
---|
445 |
|
---|
446 | Int_t Status;
|
---|
447 | Int_t M1, M2, D1, D2;
|
---|
448 |
|
---|
449 | Int_t Charge;
|
---|
450 |
|
---|
451 | Float_t Mass;
|
---|
452 |
|
---|
453 | Int_t IsPU;
|
---|
454 | Int_t IsConstituent;
|
---|
455 |
|
---|
456 | UInt_t BTag;
|
---|
457 | UInt_t TauTag;
|
---|
458 |
|
---|
459 | Float_t Eem;
|
---|
460 | Float_t Ehad;
|
---|
461 |
|
---|
462 | Float_t Edges[4];
|
---|
463 | Float_t DeltaEta;
|
---|
464 | Float_t DeltaPhi;
|
---|
465 |
|
---|
466 | TLorentzVector Momentum, Position, Area;
|
---|
467 |
|
---|
468 | Float_t Dxy;
|
---|
469 | Float_t SDxy;
|
---|
470 | Float_t Xd;
|
---|
471 | Float_t Yd;
|
---|
472 | Float_t Zd;
|
---|
473 |
|
---|
474 | // PileUpJetID variables
|
---|
475 |
|
---|
476 | Int_t NCharged;
|
---|
477 | Int_t NNeutrals;
|
---|
478 | Float_t Beta;
|
---|
479 | Float_t BetaStar;
|
---|
480 | Float_t MeanSqDeltaR;
|
---|
481 | Float_t PTD;
|
---|
482 | Float_t FracPt[5];
|
---|
483 |
|
---|
484 | // N-subjettiness variables
|
---|
485 |
|
---|
486 | Float_t Tau[5];
|
---|
487 |
|
---|
488 | static CompBase *fgCompare; //!
|
---|
489 | const CompBase *GetCompare() const { return fgCompare; }
|
---|
490 |
|
---|
491 | void AddCandidate(Candidate *object);
|
---|
492 | TObjArray *GetCandidates();
|
---|
493 |
|
---|
494 | Bool_t Overlaps(const Candidate *object) const;
|
---|
495 |
|
---|
496 | virtual void Copy(TObject &object) const;
|
---|
497 | virtual TObject *Clone(const char *newname = "") const;
|
---|
498 | virtual void Clear(Option_t* option = "");
|
---|
499 |
|
---|
500 | private:
|
---|
501 | DelphesFactory *fFactory; //!
|
---|
502 | TObjArray *fArray; //!
|
---|
503 |
|
---|
504 | void SetFactory(DelphesFactory *factory) { fFactory = factory; }
|
---|
505 |
|
---|
506 | ClassDef(Candidate, 2)
|
---|
507 | };
|
---|
508 |
|
---|
509 | #endif // DelphesClasses_h
|
---|
510 |
|
---|
511 |
|
---|