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