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