Fork me on GitHub

source: svn/trunk/classes/DelphesClasses.h@ 1177

Last change on this file since 1177 was 1123, checked in by Pavel Demin, 12 years ago

add Weight branch and Weighter module

File size: 10.1 KB
Line 
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
27class DelphesFactory;
28
29//---------------------------------------------------------------------------
30
31class Event: public TObject
32{
33public:
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
45class LHCOEvent: public Event
46{
47public:
48
49 Int_t Trigger; // trigger word
50
51 ClassDef(LHCOEvent, 1)
52};
53
54//---------------------------------------------------------------------------
55
56class LHEFEvent: public Event
57{
58public:
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
72class HepMCEvent: public Event
73{
74public:
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
101class GenParticle: public SortableObject
102{
103public:
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
110 Int_t M1; // particle 1st mother | hepevt.jmohep[number][0] - 1
111 Int_t M2; // particle 2nd mother | hepevt.jmohep[number][1] - 1
112
113 Int_t D1; // particle 1st daughter | hepevt.jdahep[number][0] - 1
114 Int_t D2; // particle last daughter | hepevt.jdahep[number][1] - 1
115
116 Int_t Charge; // particle charge
117
118 Float_t Mass; // particle mass
119
120 Float_t E; // particle energy | hepevt.phep[number][3]
121 Float_t Px; // particle momentum vector (x component) | hepevt.phep[number][0]
122 Float_t Py; // particle momentum vector (y component) | hepevt.phep[number][1]
123 Float_t Pz; // particle momentum vector (z component) | hepevt.phep[number][2]
124
125 Float_t PT; // particle transverse momentum
126 Float_t Eta; // particle pseudorapidity
127 Float_t Phi; // particle azimuthal angle
128
129 Float_t Rapidity; // particle rapidity
130
131 Float_t T; // particle vertex position (t component) | hepevt.vhep[number][3]
132 Float_t X; // particle vertex position (x component) | hepevt.vhep[number][0]
133 Float_t Y; // particle vertex position (y component) | hepevt.vhep[number][1]
134 Float_t Z; // particle vertex position (z component) | hepevt.vhep[number][2]
135
136 static CompBase *fgCompare; //!
137 const CompBase *GetCompare() const { return fgCompare; }
138
139 TLorentzVector P4();
140
141 ClassDef(GenParticle, 1)
142};
143
144//---------------------------------------------------------------------------
145
146class MissingET: public TObject
147{
148public:
149 Float_t MET; // mising transverse energy
150 Float_t Phi; // mising energy azimuthal angle
151
152 ClassDef(MissingET, 1)
153};
154
155//---------------------------------------------------------------------------
156
157class ScalarHT: public TObject
158{
159public:
160 Float_t HT; // scalar sum of transverse momenta
161
162 ClassDef(ScalarHT, 1)
163};
164
165//---------------------------------------------------------------------------
166
167class Rho: public TObject
168{
169public:
170 Float_t Rho; // rho energy density
171
172 ClassDef(Rho, 1)
173};
174
175//---------------------------------------------------------------------------
176
177class Weight: public TObject
178{
179public:
180 Float_t Weight; // weight for the event
181
182 ClassDef(Weight, 1)
183};
184
185//---------------------------------------------------------------------------
186
187class Photon: public SortableObject
188{
189public:
190
191 Float_t PT; // photon transverse momentum
192 Float_t Eta; // photon pseudorapidity
193 Float_t Phi; // photon azimuthal angle
194
195 Float_t E; // photon energy
196
197 Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
198
199 TRefArray Particles; // references to generated particles
200
201 static CompBase *fgCompare; //!
202 const CompBase *GetCompare() const { return fgCompare; }
203
204 TLorentzVector P4();
205
206 ClassDef(Photon, 2)
207};
208
209//---------------------------------------------------------------------------
210
211class Electron: public SortableObject
212{
213public:
214
215 Float_t PT; // electron transverse momentum
216 Float_t Eta; // electron pseudorapidity
217 Float_t Phi; // electron azimuthal angle
218
219 Int_t Charge; // electron charge
220
221 Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
222
223 TRef Particle; // reference to generated particle
224
225 static CompBase *fgCompare; //!
226 const CompBase *GetCompare() const { return fgCompare; }
227
228 TLorentzVector P4();
229
230 ClassDef(Electron, 2)
231};
232
233//---------------------------------------------------------------------------
234
235class Muon: public SortableObject
236{
237public:
238
239 Float_t PT; // muon transverse momentum
240 Float_t Eta; // muon pseudorapidity
241 Float_t Phi; // muon azimuthal angle
242
243 Int_t Charge; // muon charge
244
245 TRef Particle; // reference to generated particle
246
247 static CompBase *fgCompare; //!
248 const CompBase *GetCompare() const { return fgCompare; }
249
250 TLorentzVector P4();
251
252 ClassDef(Muon, 2)
253};
254
255//---------------------------------------------------------------------------
256
257class Jet: public SortableObject
258{
259public:
260
261 Float_t PT; // jet transverse momentum
262 Float_t Eta; // jet pseudorapidity
263 Float_t Phi; // jet azimuthal angle
264
265 Float_t Mass; // jet invariant mass
266
267 Float_t DeltaEta; // jet radius in pseudorapidity
268 Float_t DeltaPhi; // jet radius in azimuthal angle
269
270 UInt_t BTag; // 0 or 1 for a jet that has been tagged as containing a heavy quark
271 UInt_t TauTag; // 0 or 1 for a jet that has been tagged as a tau
272
273 Int_t Charge; // tau charge
274
275 Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
276
277 TRefArray Constituents; // references to constituents
278 TRefArray Particles; // references to generated particles
279
280 static CompBase *fgCompare; //!
281 const CompBase *GetCompare() const { return fgCompare; }
282
283 TLorentzVector P4();
284
285 ClassDef(Jet, 2)
286};
287
288//---------------------------------------------------------------------------
289
290class Track: public SortableObject
291{
292public:
293 Int_t PID; // HEP ID number
294
295 Int_t Charge; // track charge
296
297 Float_t PT; // track transverse momentum
298
299 Float_t Eta; // track pseudorapidity
300 Float_t Phi; // track azimuthal angle
301
302 Float_t EtaOuter; // track pseudorapidity at the tracker edge
303 Float_t PhiOuter; // track azimuthal angle at the tracker edge
304
305 Float_t X; // track vertex position (x component)
306 Float_t Y; // track vertex position (y component)
307 Float_t Z; // track vertex position (z component)
308
309 Float_t XOuter; // track position (x component) at the tracker edge
310 Float_t YOuter; // track position (y component) at the tracker edge
311 Float_t ZOuter; // track position (z component) at the tracker edge
312
313 TRef Particle; // reference to generated particle
314
315 static CompBase *fgCompare; //!
316 const CompBase *GetCompare() const { return fgCompare; }
317
318 TLorentzVector P4();
319
320 ClassDef(Track, 1)
321};
322
323//---------------------------------------------------------------------------
324
325class Tower: public SortableObject
326{
327public:
328 Float_t ET; // calorimeter tower transverse energy
329 Float_t Eta; // calorimeter tower pseudorapidity
330 Float_t Phi; // calorimeter tower azimuthal angle
331
332 Float_t E; // calorimeter tower energy
333
334 Float_t Eem; // calorimeter tower electromagnetic energy
335 Float_t Ehad; // calorimeter tower hadronic energy
336
337 Float_t Edges[4]; // calorimeter tower edges
338
339 TRefArray Particles; // references to generated particles
340
341 static CompBase *fgCompare; //!
342 const CompBase *GetCompare() const { return fgCompare; }
343
344 TLorentzVector P4();
345
346 ClassDef(Tower, 1)
347};
348
349//---------------------------------------------------------------------------
350
351class Candidate: public SortableObject
352{
353 friend class DelphesFactory;
354
355public:
356 Candidate();
357
358 Int_t PID;
359
360 Int_t Status;
361 Int_t M1, M2, D1, D2;
362
363 Int_t Charge;
364
365 Float_t Mass;
366
367 Int_t IsPU;
368 Int_t IsConstituent;
369
370 UInt_t BTag;
371 UInt_t TauTag;
372
373 Float_t Eem;
374 Float_t Ehad;
375
376 Float_t Edges[4];
377 Float_t DeltaEta;
378 Float_t DeltaPhi;
379
380 TLorentzVector Momentum, Position, Area;
381
382 static CompBase *fgCompare; //!
383 const CompBase *GetCompare() const { return fgCompare; }
384
385 void AddCandidate(Candidate *object);
386 TObjArray *GetCandidates();
387
388 Bool_t Overlaps(const Candidate *object) const;
389
390 virtual void Copy(TObject &object) const;
391 virtual TObject *Clone(const char *newname = "") const;
392 virtual void Clear(Option_t* option = "");
393
394private:
395 DelphesFactory *fFactory; //!
396 TObjArray *fArray; //!
397
398 void SetFactory(DelphesFactory *factory) { fFactory = factory; }
399
400 ClassDef(Candidate, 1)
401};
402
403#endif // DelphesClasses_h
404
405
Note: See TracBrowser for help on using the repository browser.