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 "TLorentzVector.h"
|
---|
35 | #include "TMatrixDSym.h"
|
---|
36 | #include "TObject.h"
|
---|
37 | #include "TRef.h"
|
---|
38 | #include "TRefArray.h"
|
---|
39 |
|
---|
40 | #include "classes/SortableObject.h"
|
---|
41 |
|
---|
42 | class DelphesFactory;
|
---|
43 |
|
---|
44 | //---------------------------------------------------------------------------
|
---|
45 |
|
---|
46 | class Event: public TObject
|
---|
47 | {
|
---|
48 | public:
|
---|
49 | Long64_t Number; // event number
|
---|
50 |
|
---|
51 | Float_t ReadTime; // read time
|
---|
52 | Float_t ProcTime; // processing time
|
---|
53 |
|
---|
54 | ClassDef(Event, 1)
|
---|
55 | };
|
---|
56 |
|
---|
57 | //---------------------------------------------------------------------------
|
---|
58 |
|
---|
59 | class LHCOEvent: public Event
|
---|
60 | {
|
---|
61 | public:
|
---|
62 | Int_t Trigger; // trigger word
|
---|
63 |
|
---|
64 | ClassDef(LHCOEvent, 1)
|
---|
65 | };
|
---|
66 |
|
---|
67 | //---------------------------------------------------------------------------
|
---|
68 |
|
---|
69 | class LHEFEvent: public Event
|
---|
70 | {
|
---|
71 | public:
|
---|
72 | Int_t ProcessID; // subprocess code for the event | hepup.IDPRUP
|
---|
73 |
|
---|
74 | Float_t Weight; // weight for the event | hepup.XWGTUP
|
---|
75 | Float_t CrossSection; // cross-section (read from init, implemented only for Wizard evgen)
|
---|
76 | Float_t ScalePDF; // scale in GeV used in the calculation of the PDFs in the event | hepup.SCALUP
|
---|
77 | Float_t AlphaQED; // value of the QED coupling used in the event | hepup.AQEDUP
|
---|
78 | Float_t AlphaQCD; // value of the QCD coupling used in the event | hepup.AQCDUP
|
---|
79 |
|
---|
80 | ClassDef(LHEFEvent, 3)
|
---|
81 | };
|
---|
82 |
|
---|
83 | //---------------------------------------------------------------------------
|
---|
84 |
|
---|
85 | class LHEFWeight: public TObject
|
---|
86 | {
|
---|
87 | public:
|
---|
88 | Int_t ID; // weight ID
|
---|
89 | Float_t Weight; // weight value
|
---|
90 |
|
---|
91 | ClassDef(LHEFWeight, 1)
|
---|
92 | };
|
---|
93 |
|
---|
94 | //---------------------------------------------------------------------------
|
---|
95 |
|
---|
96 | class HepMCEvent: public Event
|
---|
97 | {
|
---|
98 | public:
|
---|
99 | Int_t ProcessID; // unique signal process id | signal_process_id()
|
---|
100 | Int_t MPI; // number of multi parton interactions | mpi ()
|
---|
101 |
|
---|
102 | Float_t Weight; // weight for the event
|
---|
103 | Float_t CrossSection; // cross-section in pb
|
---|
104 | Float_t CrossSectionError; // cross-section error in pb
|
---|
105 |
|
---|
106 | Float_t Scale; // energy scale, see hep-ph/0109068 | event_scale()
|
---|
107 | Float_t AlphaQED; // QED coupling, see hep-ph/0109068 | alphaQED()
|
---|
108 | Float_t AlphaQCD; // QCD coupling, see hep-ph/0109068 | alphaQCD()
|
---|
109 |
|
---|
110 | Int_t ID1; // flavour code of first parton | pdf_info()->id1()
|
---|
111 | Int_t ID2; // flavour code of second parton | pdf_info()->id2()
|
---|
112 |
|
---|
113 | Float_t X1; // fraction of beam momentum carried by first parton ("beam side") | pdf_info()->x1()
|
---|
114 | Float_t X2; // fraction of beam momentum carried by second parton ("target side") | pdf_info()->x2()
|
---|
115 |
|
---|
116 | Float_t ScalePDF; // Q-scale used in evaluation of PDF's (in GeV) | pdf_info()->scalePDF()
|
---|
117 |
|
---|
118 | Float_t PDF1; // PDF (id1, x1, Q) | pdf_info()->pdf1()
|
---|
119 | Float_t PDF2; // PDF (id2, x2, Q) | pdf_info()->pdf2()
|
---|
120 |
|
---|
121 | ClassDef(HepMCEvent, 3)
|
---|
122 | };
|
---|
123 |
|
---|
124 | //---------------------------------------------------------------------------
|
---|
125 |
|
---|
126 | class GenParticle: public SortableObject
|
---|
127 | {
|
---|
128 | public:
|
---|
129 | Int_t PID; // particle HEP ID number | hepevt.idhep[number]
|
---|
130 |
|
---|
131 | Int_t Status; // particle status | hepevt.isthep[number]
|
---|
132 | Int_t IsPU; // 0 or 1 for particles from pile-up interactions
|
---|
133 |
|
---|
134 | Int_t M1; // particle 1st mother | hepevt.jmohep[number][0] - 1
|
---|
135 | Int_t M2; // particle 2nd mother | hepevt.jmohep[number][1] - 1
|
---|
136 |
|
---|
137 | Int_t D1; // particle 1st daughter | hepevt.jdahep[number][0] - 1
|
---|
138 | Int_t D2; // particle last daughter | hepevt.jdahep[number][1] - 1
|
---|
139 |
|
---|
140 | Int_t Charge; // particle charge
|
---|
141 |
|
---|
142 | Float_t Mass; // particle mass
|
---|
143 |
|
---|
144 | Float_t E; // particle energy | hepevt.phep[number][3]
|
---|
145 | Float_t Px; // particle momentum vector (x component) | hepevt.phep[number][0]
|
---|
146 | Float_t Py; // particle momentum vector (y component) | hepevt.phep[number][1]
|
---|
147 | Float_t Pz; // particle momentum vector (z component) | hepevt.phep[number][2]
|
---|
148 |
|
---|
149 | Float_t P; // particle momentum
|
---|
150 | Float_t PT; // particle transverse momentum
|
---|
151 | Float_t Eta; // particle pseudorapidity
|
---|
152 | Float_t Phi; // particle azimuthal angle
|
---|
153 |
|
---|
154 | Float_t Rapidity; // particle rapidity
|
---|
155 | Float_t CtgTheta; // particle cotangent of theta
|
---|
156 |
|
---|
157 | Float_t D0; // particle transverse impact parameter
|
---|
158 | Float_t DZ; // particle longitudinal impact parameter
|
---|
159 |
|
---|
160 | Float_t T; // particle vertex position (t component) | hepevt.vhep[number][3]
|
---|
161 | Float_t X; // particle vertex position (x component) | hepevt.vhep[number][0]
|
---|
162 | Float_t Y; // particle vertex position (y component) | hepevt.vhep[number][1]
|
---|
163 | Float_t Z; // particle vertex position (z component) | hepevt.vhep[number][2]
|
---|
164 |
|
---|
165 | static CompBase *fgCompare; //!
|
---|
166 | const CompBase *GetCompare() const { return fgCompare; }
|
---|
167 |
|
---|
168 | TLorentzVector P4() const;
|
---|
169 |
|
---|
170 | ClassDef(GenParticle, 2)
|
---|
171 | };
|
---|
172 |
|
---|
173 | //---------------------------------------------------------------------------
|
---|
174 |
|
---|
175 | class Vertex: public SortableObject
|
---|
176 | {
|
---|
177 | public:
|
---|
178 | Float_t T; // vertex position (t component)
|
---|
179 | Float_t X; // vertex position (x component)
|
---|
180 | Float_t Y; // vertex position (y component)
|
---|
181 | Float_t Z; // vertex position (z component)
|
---|
182 |
|
---|
183 | Double_t ErrorT; // vertex position error (t component)
|
---|
184 | Double_t ErrorX; // vertex position error (x component)
|
---|
185 | Double_t ErrorY; // vertex position error (y component)
|
---|
186 | Double_t ErrorZ; // vertex position error (z component)
|
---|
187 |
|
---|
188 | Int_t Index; // vertex index
|
---|
189 | Int_t NDF; // number of degrees of freedom
|
---|
190 |
|
---|
191 | Double_t Sigma; // vertex position (z component) error
|
---|
192 | Double_t SumPT2; // sum pt^2 of tracks attached to the vertex
|
---|
193 | Double_t GenSumPT2; // sum pt^2 of gen tracks attached to the vertex
|
---|
194 |
|
---|
195 | Double_t GenDeltaZ; // distance in z to closest generated vertex
|
---|
196 | Double_t BTVSumPT2; // sum pt^2 of tracks attached to the secondary vertex
|
---|
197 |
|
---|
198 | TRefArray Constituents; // references to constituents
|
---|
199 |
|
---|
200 | static CompBase *fgCompare; //!
|
---|
201 | const CompBase *GetCompare() const { return fgCompare; }
|
---|
202 |
|
---|
203 | ClassDef(Vertex, 3)
|
---|
204 | };
|
---|
205 |
|
---|
206 | //---------------------------------------------------------------------------
|
---|
207 |
|
---|
208 | class MissingET: public TObject
|
---|
209 | {
|
---|
210 | public:
|
---|
211 | Float_t MET; // mising transverse energy
|
---|
212 | Float_t Eta; // mising energy pseudorapidity
|
---|
213 | Float_t Phi; // mising energy azimuthal angle
|
---|
214 |
|
---|
215 | TLorentzVector P4() const;
|
---|
216 |
|
---|
217 | ClassDef(MissingET, 1)
|
---|
218 | };
|
---|
219 |
|
---|
220 | //---------------------------------------------------------------------------
|
---|
221 |
|
---|
222 | class ScalarHT: public TObject
|
---|
223 | {
|
---|
224 | public:
|
---|
225 | Float_t HT; // scalar sum of transverse momenta
|
---|
226 |
|
---|
227 | ClassDef(ScalarHT, 1)
|
---|
228 | };
|
---|
229 |
|
---|
230 | //---------------------------------------------------------------------------
|
---|
231 |
|
---|
232 | class Rho: public TObject
|
---|
233 | {
|
---|
234 | public:
|
---|
235 | Float_t Rho; // rho energy density
|
---|
236 | Float_t Edges[2]; // pseudorapidity range edges
|
---|
237 |
|
---|
238 | ClassDef(Rho, 1)
|
---|
239 | };
|
---|
240 |
|
---|
241 | //---------------------------------------------------------------------------
|
---|
242 |
|
---|
243 | class Weight: public TObject
|
---|
244 | {
|
---|
245 | public:
|
---|
246 | Float_t Weight; // weight for the event
|
---|
247 |
|
---|
248 | ClassDef(Weight, 1)
|
---|
249 | };
|
---|
250 |
|
---|
251 | //---------------------------------------------------------------------------
|
---|
252 |
|
---|
253 | class Photon: public SortableObject
|
---|
254 | {
|
---|
255 | public:
|
---|
256 | Float_t PT; // photon transverse momentum
|
---|
257 | Float_t Eta; // photon pseudorapidity
|
---|
258 | Float_t Phi; // photon azimuthal angle
|
---|
259 |
|
---|
260 | Float_t E; // photon energy
|
---|
261 |
|
---|
262 | Float_t T; // particle arrival time of flight
|
---|
263 |
|
---|
264 | Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
|
---|
265 |
|
---|
266 | TRefArray Particles; // references to generated particles
|
---|
267 |
|
---|
268 | Float_t IsolationVar; // isolation variable
|
---|
269 | Float_t IsolationVarRhoCorr; // isolation variable
|
---|
270 | Float_t SumPtCharged; // isolation variable
|
---|
271 | Float_t SumPtNeutral; // isolation variable
|
---|
272 | Float_t SumPtChargedPU; // isolation variable
|
---|
273 | Float_t SumPt; // isolation variable
|
---|
274 |
|
---|
275 | Int_t Status; // 1: prompt -- 2: non prompt -- 3: fake
|
---|
276 |
|
---|
277 | static CompBase *fgCompare; //!
|
---|
278 | const CompBase *GetCompare() const { return fgCompare; }
|
---|
279 |
|
---|
280 | TLorentzVector P4() const;
|
---|
281 |
|
---|
282 | ClassDef(Photon, 4)
|
---|
283 | };
|
---|
284 |
|
---|
285 | //---------------------------------------------------------------------------
|
---|
286 |
|
---|
287 | class Electron: public SortableObject
|
---|
288 | {
|
---|
289 | public:
|
---|
290 | Float_t PT; // electron transverse momentum
|
---|
291 | Float_t Eta; // electron pseudorapidity
|
---|
292 | Float_t Phi; // electron azimuthal angle
|
---|
293 |
|
---|
294 | Float_t T; // particle arrival time of flight
|
---|
295 |
|
---|
296 | Int_t Charge; // electron charge
|
---|
297 |
|
---|
298 | Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
|
---|
299 |
|
---|
300 | TRef Particle; // reference to generated particle
|
---|
301 |
|
---|
302 | Float_t IsolationVar; // isolation variable
|
---|
303 | Float_t IsolationVarRhoCorr; // isolation variable
|
---|
304 | Float_t SumPtCharged; // isolation variable
|
---|
305 | Float_t SumPtNeutral; // isolation variable
|
---|
306 | Float_t SumPtChargedPU; // isolation variable
|
---|
307 | Float_t SumPt; // isolation variable
|
---|
308 |
|
---|
309 | Float_t D0; // track transverse impact parameter
|
---|
310 | Float_t DZ; // track longitudinal impact parameter
|
---|
311 | Float_t ErrorD0; // track transverse impact parameter error
|
---|
312 | Float_t ErrorDZ; // track longitudinal impact parameter error
|
---|
313 |
|
---|
314 | static CompBase *fgCompare; //!
|
---|
315 | const CompBase *GetCompare() const { return fgCompare; }
|
---|
316 |
|
---|
317 | TLorentzVector P4() const;
|
---|
318 |
|
---|
319 | ClassDef(Electron, 4)
|
---|
320 | };
|
---|
321 |
|
---|
322 | //---------------------------------------------------------------------------
|
---|
323 |
|
---|
324 | class Muon: public SortableObject
|
---|
325 | {
|
---|
326 | public:
|
---|
327 | Float_t PT; // muon transverse momentum
|
---|
328 | Float_t Eta; // muon pseudorapidity
|
---|
329 | Float_t Phi; // muon azimuthal angle
|
---|
330 |
|
---|
331 | Float_t T; // particle arrival time of flight
|
---|
332 |
|
---|
333 | Int_t Charge; // muon charge
|
---|
334 |
|
---|
335 | TRef Particle; // reference to generated particle
|
---|
336 |
|
---|
337 | Float_t IsolationVar; // isolation variable
|
---|
338 | Float_t IsolationVarRhoCorr; // isolation variable
|
---|
339 | Float_t SumPtCharged; // isolation variable
|
---|
340 | Float_t SumPtNeutral; // isolation variable
|
---|
341 | Float_t SumPtChargedPU; // isolation variable
|
---|
342 | Float_t SumPt; // isolation variable
|
---|
343 |
|
---|
344 | Float_t D0; // track transverse impact parameter
|
---|
345 | Float_t DZ; // track longitudinal impact parameter
|
---|
346 | Float_t ErrorD0; // track transverse impact parameter error
|
---|
347 | Float_t ErrorDZ; // track longitudinal impact parameter error
|
---|
348 |
|
---|
349 | static CompBase *fgCompare; //!
|
---|
350 | const CompBase *GetCompare() const { return fgCompare; }
|
---|
351 |
|
---|
352 | TLorentzVector P4() const;
|
---|
353 |
|
---|
354 | ClassDef(Muon, 4)
|
---|
355 | };
|
---|
356 |
|
---|
357 | //---------------------------------------------------------------------------
|
---|
358 |
|
---|
359 | class Jet: public SortableObject
|
---|
360 | {
|
---|
361 | public:
|
---|
362 | Float_t PT; // jet transverse momentum
|
---|
363 | Float_t Eta; // jet pseudorapidity
|
---|
364 | Float_t Phi; // jet azimuthal angle
|
---|
365 |
|
---|
366 | Float_t T; //particle arrival time of flight
|
---|
367 |
|
---|
368 | Float_t Mass; // jet invariant mass
|
---|
369 |
|
---|
370 | Float_t DeltaEta; // jet radius in pseudorapidity
|
---|
371 | Float_t DeltaPhi; // jet radius in azimuthal angle
|
---|
372 |
|
---|
373 | UInt_t Flavor; // jet flavor
|
---|
374 | UInt_t FlavorAlgo; // jet flavor
|
---|
375 | UInt_t FlavorPhys; // jet flavor
|
---|
376 |
|
---|
377 | UInt_t BTag; // 0 or 1 for a jet that has been tagged as containing a heavy quark
|
---|
378 | UInt_t BTagAlgo; // 0 or 1 for a jet that has been tagged as containing a heavy quark
|
---|
379 | UInt_t BTagPhys; // 0 or 1 for a jet that has been tagged as containing a heavy quark
|
---|
380 |
|
---|
381 | UInt_t TauTag; // 0 or 1 for a jet that has been tagged as a tau
|
---|
382 | Float_t TauWeight; // probability for jet to be identified as tau
|
---|
383 |
|
---|
384 | Int_t Charge; // tau charge
|
---|
385 |
|
---|
386 | Float_t EhadOverEem; // ratio of the hadronic versus electromagnetic energy deposited in the calorimeter
|
---|
387 |
|
---|
388 | Int_t NCharged; // number of charged constituents
|
---|
389 | Int_t NNeutrals; // number of neutral constituents
|
---|
390 |
|
---|
391 | Float_t NeutralEnergyFraction; // charged energy fraction
|
---|
392 | Float_t ChargedEnergyFraction; // neutral energy fraction
|
---|
393 |
|
---|
394 | Float_t Beta; // (sum pt of charged pile-up constituents)/(sum pt of charged constituents)
|
---|
395 | Float_t BetaStar; // (sum pt of charged constituents coming from hard interaction)/(sum pt of charged constituents)
|
---|
396 | Float_t MeanSqDeltaR; // average distance (squared) between constituent and jet weighted by pt (squared) of constituent
|
---|
397 | Float_t PTD; // average pt between constituent and jet weighted by pt of constituent
|
---|
398 | Float_t FracPt[5]; // (sum pt of constituents within a ring 0.1*i < DeltaR < 0.1*(i+1))/(sum pt of constituents)
|
---|
399 |
|
---|
400 | Float_t Tau[5]; // N-subjettiness
|
---|
401 |
|
---|
402 | TLorentzVector SoftDroppedJet;
|
---|
403 | TLorentzVector SoftDroppedSubJet1;
|
---|
404 | TLorentzVector SoftDroppedSubJet2;
|
---|
405 |
|
---|
406 | TLorentzVector TrimmedP4[5]; // first entry (i = 0) is the total Trimmed Jet 4-momenta and from i = 1 to 4 are the trimmed subjets 4-momenta
|
---|
407 | TLorentzVector PrunedP4[5]; // first entry (i = 0) is the total Pruned Jet 4-momenta and from i = 1 to 4 are the pruned subjets 4-momenta
|
---|
408 | TLorentzVector SoftDroppedP4[5]; // first entry (i = 0) is the total SoftDropped Jet 4-momenta and from i = 1 to 4 are the pruned subjets 4-momenta
|
---|
409 |
|
---|
410 | Int_t NSubJetsTrimmed; // number of subjets trimmed
|
---|
411 | Int_t NSubJetsPruned; // number of subjets pruned
|
---|
412 | Int_t NSubJetsSoftDropped; // number of subjets soft-dropped
|
---|
413 |
|
---|
414 | Double_t ExclYmerge23;
|
---|
415 | Double_t ExclYmerge34;
|
---|
416 | Double_t ExclYmerge45;
|
---|
417 | Double_t ExclYmerge56;
|
---|
418 |
|
---|
419 | TRefArray Constituents; // references to constituents
|
---|
420 | TRefArray Particles; // references to generated particles
|
---|
421 |
|
---|
422 | static CompBase *fgCompare; //!
|
---|
423 | const CompBase *GetCompare() const { return fgCompare; }
|
---|
424 |
|
---|
425 | TLorentzVector P4() const;
|
---|
426 | TLorentzVector Area;
|
---|
427 |
|
---|
428 | ClassDef(Jet, 4)
|
---|
429 | };
|
---|
430 |
|
---|
431 | //---------------------------------------------------------------------------
|
---|
432 |
|
---|
433 | class Track: public SortableObject
|
---|
434 | {
|
---|
435 | public:
|
---|
436 | Int_t PID; // HEP ID number
|
---|
437 |
|
---|
438 | Int_t Charge; // track charge
|
---|
439 |
|
---|
440 | Float_t P; // track momentum
|
---|
441 | Float_t PT; // track transverse momentum
|
---|
442 | Float_t Eta; // track pseudorapidity
|
---|
443 | Float_t Phi; // track azimuthal angle
|
---|
444 | Float_t CtgTheta; // track cotangent of theta
|
---|
445 |
|
---|
446 | Float_t EtaOuter; // track pseudorapidity at the tracker edge
|
---|
447 | Float_t PhiOuter; // track azimuthal angle at the tracker edge
|
---|
448 |
|
---|
449 | Float_t T; // track vertex position (t component)
|
---|
450 | Float_t X; // track vertex position (x component)
|
---|
451 | Float_t Y; // track vertex position (y component)
|
---|
452 | Float_t Z; // track vertex position (z component)
|
---|
453 |
|
---|
454 | Float_t TOuter; // track position (t component) at the tracker edge
|
---|
455 | Float_t XOuter; // track position (x component) at the tracker edge
|
---|
456 | Float_t YOuter; // track position (y component) at the tracker edge
|
---|
457 | Float_t ZOuter; // track position (z component) at the tracker edge
|
---|
458 |
|
---|
459 | Float_t Xd; // X coordinate of point of closest approach to vertex
|
---|
460 | Float_t Yd; // Y coordinate of point of closest approach to vertex
|
---|
461 | Float_t Zd; // Z coordinate of point of closest approach to vertex
|
---|
462 |
|
---|
463 | Float_t L; // track path length
|
---|
464 | Float_t D0; // track transverse impact parameter
|
---|
465 | Float_t DZ; // track longitudinal impact parameter
|
---|
466 |
|
---|
467 | Float_t ErrorP; // track momentum error
|
---|
468 | Float_t ErrorPT; // track transverse momentum error
|
---|
469 | Float_t ErrorPhi; // track azimuthal angle error
|
---|
470 | Float_t ErrorCtgTheta; // track cotangent of theta error
|
---|
471 |
|
---|
472 | Float_t ErrorT; // time measurement error
|
---|
473 | Float_t ErrorD0; // track transverse impact parameter error
|
---|
474 | Float_t ErrorDZ; // track longitudinal impact parameter error
|
---|
475 | Float_t ErrorC; // track curvature error
|
---|
476 |
|
---|
477 | // track covariance off-diagonal terms
|
---|
478 | Float_t ErrorD0Phi;
|
---|
479 | Float_t ErrorD0C;
|
---|
480 | Float_t ErrorD0DZ;
|
---|
481 | Float_t ErrorD0CtgTheta;
|
---|
482 | Float_t ErrorPhiC;
|
---|
483 | Float_t ErrorPhiDZ;
|
---|
484 | Float_t ErrorPhiCtgTheta ;
|
---|
485 | Float_t ErrorCDZ;
|
---|
486 | Float_t ErrorCCtgTheta;
|
---|
487 | Float_t ErrorDZCtgTheta;
|
---|
488 |
|
---|
489 | TRef Particle; // reference to generated particle
|
---|
490 |
|
---|
491 | Int_t VertexIndex; // reference to vertex
|
---|
492 |
|
---|
493 | static CompBase *fgCompare; //!
|
---|
494 | const CompBase *GetCompare() const { return fgCompare; }
|
---|
495 |
|
---|
496 | TLorentzVector P4() const;
|
---|
497 | TMatrixDSym CovarianceMatrix() const;
|
---|
498 |
|
---|
499 | ClassDef(Track, 3)
|
---|
500 | };
|
---|
501 |
|
---|
502 | //---------------------------------------------------------------------------
|
---|
503 |
|
---|
504 | class Tower: public SortableObject
|
---|
505 | {
|
---|
506 | public:
|
---|
507 | Float_t ET; // calorimeter tower transverse energy
|
---|
508 | Float_t Eta; // calorimeter tower pseudorapidity
|
---|
509 | Float_t Phi; // calorimeter tower azimuthal angle
|
---|
510 |
|
---|
511 | Float_t E; // calorimeter tower energy
|
---|
512 |
|
---|
513 | Float_t T; // ecal deposit time, averaged by sqrt(EM energy) over all particles, not smeared
|
---|
514 | Int_t NTimeHits; // number of hits contributing to time measurement
|
---|
515 |
|
---|
516 | Float_t Eem; // calorimeter tower electromagnetic energy
|
---|
517 | Float_t Ehad; // calorimeter tower hadronic energy
|
---|
518 |
|
---|
519 | Float_t Edges[4]; // calorimeter tower edges
|
---|
520 |
|
---|
521 | TRefArray Particles; // references to generated particles
|
---|
522 |
|
---|
523 | static CompBase *fgCompare; //!
|
---|
524 | const CompBase *GetCompare() const { return fgCompare; }
|
---|
525 |
|
---|
526 | TLorentzVector P4() const;
|
---|
527 |
|
---|
528 | ClassDef(Tower, 2)
|
---|
529 | };
|
---|
530 |
|
---|
531 | //---------------------------------------------------------------------------
|
---|
532 |
|
---|
533 | class ParticleFlowCandidate: public SortableObject
|
---|
534 | {
|
---|
535 |
|
---|
536 | public:
|
---|
537 | Int_t PID; // HEP ID number
|
---|
538 |
|
---|
539 | Int_t Charge; // track charge
|
---|
540 |
|
---|
541 | Float_t E; // reconstructed energy [GeV]
|
---|
542 | Float_t P; // track momentum
|
---|
543 | Float_t PT; // track transverse momentum
|
---|
544 | Float_t Eta; // track pseudorapidity
|
---|
545 | Float_t Phi; // track azimuthal angle
|
---|
546 | Float_t CtgTheta; // track cotangent of theta
|
---|
547 |
|
---|
548 | Float_t EtaOuter; // track pseudorapidity at the tracker edge
|
---|
549 | Float_t PhiOuter; // track azimuthal angle at the tracker edge
|
---|
550 |
|
---|
551 | Float_t T; // track vertex position (t component)
|
---|
552 | Float_t X; // track vertex position (x component)
|
---|
553 | Float_t Y; // track vertex position (y component)
|
---|
554 | Float_t Z; // track vertex position (z component)
|
---|
555 |
|
---|
556 | Float_t TOuter; // track position (t component) at the tracker edge
|
---|
557 | Float_t XOuter; // track position (x component) at the tracker edge
|
---|
558 | Float_t YOuter; // track position (y component) at the tracker edge
|
---|
559 | Float_t ZOuter; // track position (z component) at the tracker edge
|
---|
560 |
|
---|
561 | Float_t Xd; // X coordinate of point of closest approach to vertex
|
---|
562 | Float_t Yd; // Y coordinate of point of closest approach to vertex
|
---|
563 | Float_t Zd; // Z coordinate of point of closest approach to vertex
|
---|
564 |
|
---|
565 | Float_t L; // track path length
|
---|
566 | Float_t D0; // track transverse impact parameter
|
---|
567 | Float_t DZ; // track longitudinal impact parameter
|
---|
568 |
|
---|
569 | Float_t ErrorP; // track momentum error
|
---|
570 | Float_t ErrorPT; // track transverse momentum error
|
---|
571 | Float_t ErrorPhi; // track azimuthal angle error
|
---|
572 | Float_t ErrorCtgTheta; // track cotangent of theta error
|
---|
573 |
|
---|
574 | Float_t ErrorT; // time measurement error
|
---|
575 | Float_t ErrorD0; // track transverse impact parameter error
|
---|
576 | Float_t ErrorDZ; // track longitudinal impact parameter error
|
---|
577 | Float_t ErrorC; // track curvature error
|
---|
578 |
|
---|
579 | // track covariance off-diagonal terms
|
---|
580 | Float_t ErrorD0Phi;
|
---|
581 | Float_t ErrorD0C;
|
---|
582 | Float_t ErrorD0DZ;
|
---|
583 | Float_t ErrorD0CtgTheta;
|
---|
584 | Float_t ErrorPhiC;
|
---|
585 | Float_t ErrorPhiDZ;
|
---|
586 | Float_t ErrorPhiCtgTheta ;
|
---|
587 | Float_t ErrorCDZ;
|
---|
588 | Float_t ErrorCCtgTheta;
|
---|
589 | Float_t ErrorDZCtgTheta;
|
---|
590 |
|
---|
591 | Int_t VertexIndex; // reference to vertex
|
---|
592 |
|
---|
593 | static CompBase *fgCompare; //!
|
---|
594 | const CompBase *GetCompare() const { return fgCompare; }
|
---|
595 |
|
---|
596 | TLorentzVector P4() const;
|
---|
597 | TMatrixDSym CovarianceMatrix() const;
|
---|
598 |
|
---|
599 | Int_t NTimeHits; // number of hits contributing to time measurement
|
---|
600 |
|
---|
601 | Float_t Eem; // calorimeter tower electromagnetic energy
|
---|
602 | Float_t Ehad; // calorimeter tower hadronic energy
|
---|
603 |
|
---|
604 | Float_t Edges[4]; // calorimeter tower edges
|
---|
605 |
|
---|
606 | TRefArray Particles; // references to generated particles
|
---|
607 |
|
---|
608 | ClassDef(ParticleFlowCandidate, 2)
|
---|
609 |
|
---|
610 | };
|
---|
611 |
|
---|
612 | //---------------------------------------------------------------------------
|
---|
613 |
|
---|
614 | class HectorHit: public SortableObject
|
---|
615 | {
|
---|
616 | public:
|
---|
617 | Float_t E; // reconstructed energy [GeV]
|
---|
618 |
|
---|
619 | Float_t Tx; // angle of the momentum in the horizontal (x,z) plane [urad]
|
---|
620 | Float_t Ty; // angle of the momentum in the verical (y,z) plane [urad]
|
---|
621 |
|
---|
622 | Float_t T; // time of flight to the detector [s]
|
---|
623 |
|
---|
624 | Float_t X; // horizontal distance to the beam [um]
|
---|
625 | Float_t Y; // vertical distance to the beam [um]
|
---|
626 | Float_t S; // distance to the interaction point [m]
|
---|
627 |
|
---|
628 | TRef Particle; // reference to generated particle
|
---|
629 |
|
---|
630 | static CompBase *fgCompare; //!
|
---|
631 | const CompBase *GetCompare() const { return fgCompare; }
|
---|
632 |
|
---|
633 | ClassDef(HectorHit, 1)
|
---|
634 | };
|
---|
635 |
|
---|
636 | //---------------------------------------------------------------------------
|
---|
637 |
|
---|
638 | class Candidate: public SortableObject
|
---|
639 | {
|
---|
640 | friend class DelphesFactory;
|
---|
641 |
|
---|
642 | public:
|
---|
643 | Candidate();
|
---|
644 |
|
---|
645 | Int_t PID;
|
---|
646 |
|
---|
647 | Int_t Status;
|
---|
648 | Int_t M1, M2, D1, D2;
|
---|
649 |
|
---|
650 | Int_t Charge;
|
---|
651 |
|
---|
652 | Float_t Mass;
|
---|
653 |
|
---|
654 | Int_t IsPU;
|
---|
655 | Int_t IsRecoPU;
|
---|
656 |
|
---|
657 | Int_t IsConstituent;
|
---|
658 | Int_t IsFromConversion;
|
---|
659 |
|
---|
660 | UInt_t Flavor;
|
---|
661 | UInt_t FlavorAlgo;
|
---|
662 | UInt_t FlavorPhys;
|
---|
663 |
|
---|
664 | UInt_t BTag;
|
---|
665 | UInt_t BTagAlgo;
|
---|
666 | UInt_t BTagPhys;
|
---|
667 |
|
---|
668 | UInt_t TauTag;
|
---|
669 | Float_t TauWeight;
|
---|
670 |
|
---|
671 | Float_t Eem;
|
---|
672 | Float_t Ehad;
|
---|
673 |
|
---|
674 | Float_t Edges[4];
|
---|
675 | Float_t DeltaEta;
|
---|
676 | Float_t DeltaPhi;
|
---|
677 |
|
---|
678 | TLorentzVector Momentum, Position, InitialPosition, PositionError, Area;
|
---|
679 |
|
---|
680 | Float_t L; // path length
|
---|
681 | Float_t DZ;
|
---|
682 | Float_t ErrorDZ;
|
---|
683 | Float_t ErrorT; // path length
|
---|
684 | Float_t D0;
|
---|
685 | Float_t ErrorD0;
|
---|
686 | Float_t C;
|
---|
687 | Float_t ErrorC;
|
---|
688 | Float_t P;
|
---|
689 | Float_t ErrorP;
|
---|
690 | Float_t PT;
|
---|
691 | Float_t ErrorPT;
|
---|
692 | Float_t CtgTheta;
|
---|
693 | Float_t ErrorCtgTheta;
|
---|
694 | Float_t Phi;
|
---|
695 | Float_t ErrorPhi;
|
---|
696 |
|
---|
697 | Float_t Xd;
|
---|
698 | Float_t Yd;
|
---|
699 | Float_t Zd;
|
---|
700 |
|
---|
701 | // tracking resolution
|
---|
702 |
|
---|
703 | Float_t TrackResolution;
|
---|
704 |
|
---|
705 | // PileUpJetID variables
|
---|
706 |
|
---|
707 | Int_t NCharged;
|
---|
708 | Int_t NNeutrals;
|
---|
709 | Float_t Beta;
|
---|
710 | Float_t BetaStar;
|
---|
711 | Float_t MeanSqDeltaR;
|
---|
712 | Float_t PTD;
|
---|
713 | Float_t FracPt[5];
|
---|
714 | Float_t NeutralEnergyFraction; // charged energy fraction
|
---|
715 | Float_t ChargedEnergyFraction; // neutral energy fraction
|
---|
716 |
|
---|
717 |
|
---|
718 | // Timing information
|
---|
719 |
|
---|
720 | Int_t NTimeHits;
|
---|
721 | std::vector<std::pair<Float_t, Float_t> > ECalEnergyTimePairs;
|
---|
722 |
|
---|
723 | // Isolation variables
|
---|
724 |
|
---|
725 | Float_t IsolationVar;
|
---|
726 | Float_t IsolationVarRhoCorr;
|
---|
727 | Float_t SumPtCharged;
|
---|
728 | Float_t SumPtNeutral;
|
---|
729 | Float_t SumPtChargedPU;
|
---|
730 | Float_t SumPt;
|
---|
731 |
|
---|
732 | // ACTS compliant 6x6 track covariance (D0, phi, Curvature, dz, ctg(theta))
|
---|
733 |
|
---|
734 | TMatrixDSym TrackCovariance;
|
---|
735 |
|
---|
736 | // vertex variables
|
---|
737 |
|
---|
738 | Int_t ClusterIndex;
|
---|
739 | Int_t ClusterNDF;
|
---|
740 | Double_t ClusterSigma;
|
---|
741 | Double_t SumPT2;
|
---|
742 | Double_t BTVSumPT2;
|
---|
743 | Double_t GenDeltaZ;
|
---|
744 | Double_t GenSumPT2;
|
---|
745 |
|
---|
746 | // N-subjettiness variables
|
---|
747 |
|
---|
748 | Float_t Tau[5];
|
---|
749 |
|
---|
750 | // Other Substructure variables
|
---|
751 |
|
---|
752 | TLorentzVector SoftDroppedJet;
|
---|
753 | TLorentzVector SoftDroppedSubJet1;
|
---|
754 | TLorentzVector SoftDroppedSubJet2;
|
---|
755 |
|
---|
756 | TLorentzVector TrimmedP4[5]; // first entry (i = 0) is the total Trimmed Jet 4-momenta and from i = 1 to 4 are the trimmed subjets 4-momenta
|
---|
757 | TLorentzVector PrunedP4[5]; // first entry (i = 0) is the total Pruned Jet 4-momenta and from i = 1 to 4 are the pruned subjets 4-momenta
|
---|
758 | TLorentzVector SoftDroppedP4[5]; // first entry (i = 0) is the total SoftDropped Jet 4-momenta and from i = 1 to 4 are the pruned subjets 4-momenta
|
---|
759 |
|
---|
760 | Int_t NSubJetsTrimmed; // number of subjets trimmed
|
---|
761 | Int_t NSubJetsPruned; // number of subjets pruned
|
---|
762 | Int_t NSubJetsSoftDropped; // number of subjets soft-dropped
|
---|
763 |
|
---|
764 | // Exclusive clustering variables
|
---|
765 | Double_t ExclYmerge23;
|
---|
766 | Double_t ExclYmerge34;
|
---|
767 | Double_t ExclYmerge45;
|
---|
768 | Double_t ExclYmerge56;
|
---|
769 |
|
---|
770 | // event characteristics variables
|
---|
771 | Double_t ParticleDensity; // particle multiplicity density in the proximity of the particle
|
---|
772 |
|
---|
773 | static CompBase *fgCompare; //!
|
---|
774 | const CompBase *GetCompare() const { return fgCompare; }
|
---|
775 |
|
---|
776 | void AddCandidate(Candidate *object);
|
---|
777 | TObjArray *GetCandidates();
|
---|
778 |
|
---|
779 | Bool_t Overlaps(const Candidate *object) const;
|
---|
780 |
|
---|
781 | virtual void Copy(TObject &object) const;
|
---|
782 | virtual TObject *Clone(const char *newname = "") const;
|
---|
783 | virtual void Clear(Option_t *option = "");
|
---|
784 |
|
---|
785 | private:
|
---|
786 | DelphesFactory *fFactory; //!
|
---|
787 | TObjArray *fArray; //!
|
---|
788 |
|
---|
789 | void SetFactory(DelphesFactory *factory) { fFactory = factory; }
|
---|
790 |
|
---|
791 | ClassDef(Candidate, 6)
|
---|
792 | };
|
---|
793 |
|
---|
794 | #endif // DelphesClasses_h
|
---|