eic-smear  1.0.3
A collection of ROOT classes for Monte Carlo events and a fast-smearing code simulating detector effects for the Electron-Ion Collider task force
ParticleMC.h
Go to the documentation of this file.
1 
10 #ifndef INCLUDE_EICSMEAR_ERHIC_PARTICLEMC_H_
11 #define INCLUDE_EICSMEAR_ERHIC_PARTICLEMC_H_
12 
13 #include <string>
14 
15 #include <TLorentzVector.h>
16 #include <TRef.h>
17 #include <TVector3.h>
18 
19 #include "eicsmear/erhic/Pid.h"
21 
22 namespace erhic {
23 
24 class EventMC;
28 class ParticleMC : public VirtualParticle {
29  public:
36  explicit ParticleMC(const std::string& = "");
37 
41  virtual ~ParticleMC();
42 
49  virtual void Print(Option_t* = "") const;
50 
54  virtual UInt_t GetIndex() const;
55 
61  virtual UShort_t GetStatus() const;
62 
66  virtual UShort_t GetParentIndex() const;
67 
75  virtual const ParticleMC* GetParent() const;
76 
81  virtual UShort_t GetChild1Index() const;
82 
87  virtual UShort_t GetChildNIndex() const;
88 
93  virtual UInt_t GetNChildren() const;
94 
105  virtual const ParticleMC* GetChild(UShort_t) const;
106 
113  virtual Bool_t HasChild(Int_t) const;
114 
118  virtual Double_t GetPx() const;
119 
123  virtual Double_t GetPy() const;
124 
128  virtual Double_t GetPz() const;
129 
133  virtual Double_t GetM() const;
134 
138  virtual Double_t GetPt() const;
139 
144  virtual TVector3 GetVertex() const;
145 
150  virtual Pid GetParentId() const;
151 
155  virtual Double_t GetP() const;
156 
160  virtual Double_t GetTheta() const;
161 
165  virtual Double_t GetPhi() const;
166 
170  virtual Double_t GetRapidity() const;
171 
175  virtual Double_t GetEta() const;
176 
181  virtual Double_t GetZ() const;
182 
187  virtual Double_t GetXFeynman() const;
188 
194  virtual Double_t GetThetaVsGamma() const;
195 
200  virtual Double_t GetPtVsGamma() const;
201 
207  const EventMC* GetEvent() const;
208 
212  void SetEvent(EventMC* event);
213 
217  virtual TLorentzVector Get4Vector() const;
218 
222  virtual TLorentzVector PxPyPzE() const { return Get4Vector(); }
223 
241  virtual TLorentzVector Get4VectorInHadronBosonFrame() const;
242 
246  virtual Double_t GetE() const;
247 
248  virtual void SetE(Double_t);
249 
250  virtual void SetM(Double_t);
251 
252  virtual void SetP(Double_t);
253 
254  virtual void SetPt(Double_t);
255 
256  virtual void SetPz(Double_t);
257 
258  virtual void SetPhi(Double_t);
259 
260  virtual void SetTheta(Double_t);
261 
262  virtual void SetStatus(UShort_t);
263 
267  virtual Pid Id() const;
268 
272  virtual Pid GetPdgCode() const { return Id(); }
273 
287  virtual void ComputeDerivedQuantities();
288 
304 
310  virtual void SetIndex(int i) { I = i; }
311 
314  virtual void SetStatus(int i) { KS = i; }
315 
319  virtual void SetId(int i) { id = i; }
320 
323  virtual void SetParentIndex(int i) { orig = i; }
324 
327  virtual void SetChild1Index(int i) { daughter = i; }
328 
331  virtual void SetChildNIndex(int i) { ldaughter = i; }
332 
337  virtual void Set4Vector(const TLorentzVector&);
338 
342  virtual void SetVertex(const TVector3&);
343 
347  virtual void SetParentId(int i) { parentId = i; }
348 
349 // protected:
350 
351  UShort_t I;
352  UShort_t KS;
353  Int_t id;
354  UShort_t orig;
355  UShort_t daughter;
356  UShort_t ldaughter;
357 
358  Double32_t px;
359  Double32_t py;
360  Double32_t pz;
361 
362  // Reducing file size.
363  // The following variables could be removed, at the cost of increased
364  // running time to calculate them on-the-fly:
365  // p - compute from px/y/z
366  // m - accessible via PID with TDatabasePDG
367  // E - compute from p and m
368  // parentId - accessible via event record and parent index IF a
369  // persistent pointer to the containing event can be stored
370  // daughter - IF daughter particles always start immediately after this one
371  // pt - from px, py
372  // theta - from pt, pz
373  // phi - from px, py
374  // rapidity - from E, pz
375  // eta - from theta
376  // xFeynman - if event is accessible
377  // thetaGamma - duplicated in pHadronBoson
378  // ptVsGamma - duplicated in pHadronBoson
379  // phiPrf - duplicated in pHadronBoson
380 
381  Double32_t E;
382  Double32_t m;
383  Double32_t pt;
384  Double32_t xv;
385  Double32_t yv;
386  Double32_t zv;
387 
388  // Can be deprecated if parent particle itself is available
389  Int_t parentId;
390 
391  Double32_t p;
392  Double32_t theta;
393  Double32_t phi;
394  Double32_t rapidity;
395  Double32_t eta;
396  Double32_t z;
397  Double32_t xFeynman;
399  Double32_t thetaGamma;
400  Double32_t ptVsGamma;
402  Double32_t phiPrf;
404 
406  TRef event;
407 
409  ClassDef(ParticleMC, 1)
410 };
411 
412 inline UInt_t ParticleMC::GetIndex() const {
413  return I;
414 }
415 
416 inline UShort_t ParticleMC::GetStatus() const {
417  return KS;
418 }
419 
420 inline UShort_t ParticleMC::GetParentIndex() const {
421  return orig;
422 }
423 
424 inline UShort_t ParticleMC::GetChild1Index() const {
425  return daughter;
426 }
427 
428 inline UShort_t ParticleMC::GetChildNIndex() const {
429  return ldaughter;
430 }
431 
432 inline Double_t ParticleMC::GetPx() const {
433  return px;
434 }
435 
436 inline Double_t ParticleMC::GetPy() const {
437  return py;
438 }
439 
440 inline Double_t ParticleMC::GetPz() const {
441  return pz;
442 }
443 
444 inline Double_t ParticleMC::GetM() const {
445  return m;
446 }
447 
448 inline Double_t ParticleMC::GetPt() const {
449  return pt;
450 }
451 
452 inline TVector3 ParticleMC::GetVertex() const {
453  return TVector3(xv, yv, zv);
454 }
455 
457  return erhic::Pid(parentId);
458 }
459 
460 inline Double_t ParticleMC::GetP() const {
461  return p;
462 }
463 
464 inline Double_t ParticleMC::GetTheta() const {
465  return theta;
466 }
467 
468 inline Double_t ParticleMC::GetPhi() const {
469  return phi;
470 }
471 
472 inline Double_t ParticleMC::GetRapidity() const {
473  return rapidity;
474 }
475 
476 inline Double_t ParticleMC::GetEta() const {
477  return eta;
478 }
479 
480 inline Double_t ParticleMC::GetZ() const {
481  return z;
482 }
483 
484 inline Double_t ParticleMC::GetXFeynman() const {
485  return xFeynman;
486 }
487 
488 inline Double_t ParticleMC::GetThetaVsGamma() const {
489  return thetaGamma;
490 }
491 
492 inline Double_t ParticleMC::GetPtVsGamma() const {
493  return ptVsGamma;
494 }
495 
496 inline Pid ParticleMC::Id() const {
497  return Pid(id);
498 }
499 
500 inline UInt_t ParticleMC::GetNChildren() const {
501  if (0 == daughter) return 0;
502  if (0 == ldaughter) return 1;
503  return ldaughter - daughter + 1;
504 }
505 
506 inline Double_t ParticleMC::GetE() const {
507  return E;
508 }
509 
510 inline void ParticleMC::SetE(Double_t e) {
511  E = e;
512 }
513 
514 inline void ParticleMC::SetM(Double_t mass) {
515  m = mass;
516 }
517 
518 inline void ParticleMC::SetP(Double_t momentum) {
519  p = momentum;
520 }
521 
522 inline void ParticleMC::SetPt(Double_t momentum) {
523  pt = momentum;
524 }
525 
526 inline void ParticleMC::SetPz(Double_t momentum) {
527  pz = momentum;
528 }
529 
530 inline void ParticleMC::SetPhi(Double_t value) {
531  phi = value;
532 }
533 
534 inline void ParticleMC::SetTheta(Double_t value) {
535  theta = value;
536 }
537 
538 inline void ParticleMC::SetStatus(UShort_t status) {
539  KS = status;
540 }
541 
542 } // namespace erhic
543 
544 #endif // INCLUDE_EICSMEAR_ERHIC_PARTICLEMC_H_
Double32_t thetaGamma
Definition: ParticleMC.h:399
Double32_t m
Invariant mass of particle.
Definition: ParticleMC.h:382
UShort_t ldaughter
I of last child particle.
Definition: ParticleMC.h:356
virtual TVector3 GetVertex() const
Definition: ParticleMC.h:452
virtual Double_t GetThetaVsGamma() const
Definition: ParticleMC.h:488
virtual Double_t GetRapidity() const
Definition: ParticleMC.h:472
Double32_t px
x component of particle momentum
Definition: ParticleMC.h:358
Double32_t zv
z coordinate of particle production vertex
Definition: ParticleMC.h:386
virtual Double_t GetM() const
Definition: ParticleMC.h:444
virtual UShort_t GetChild1Index() const
Definition: ParticleMC.h:424
Double32_t E
Energy of particle.
Definition: ParticleMC.h:381
void SetEvent(EventMC *event)
Definition: ParticleMC.cxx:280
virtual void SetChildNIndex(int i)
Definition: ParticleMC.h:331
Double32_t eta
Pseudorapidity of particle.
Definition: ParticleMC.h:395
Double32_t ptVsGamma
Definition: ParticleMC.h:401
Double32_t py
y component of particle momentum
Definition: ParticleMC.h:359
virtual Double_t GetZ() const
Definition: ParticleMC.h:480
UShort_t I
Particle index in event.
Definition: ParticleMC.h:351
virtual Double_t GetPz() const
Definition: ParticleMC.h:440
virtual void ComputeDerivedQuantities()
Definition: ParticleMC.cxx:110
virtual Pid GetParentId() const
Definition: ParticleMC.h:456
Double32_t xv
x coordinate of particle production vertex
Definition: ParticleMC.h:384
virtual TLorentzVector Get4Vector() const
Definition: ParticleMC.cxx:176
ParticleMC(const std::string &="")
Definition: ParticleMC.cxx:55
virtual Double_t GetPtVsGamma() const
Definition: ParticleMC.h:492
const EventMC * GetEvent() const
Definition: ParticleMC.cxx:180
virtual void SetId(int i)
Definition: ParticleMC.h:319
virtual void Print(Option_t *="") const
Definition: ParticleMC.cxx:103
UShort_t orig
I of parent particle.
Definition: ParticleMC.h:354
Double32_t pz
z component of particle momentum
Definition: ParticleMC.h:360
virtual UShort_t GetStatus() const
Definition: ParticleMC.h:416
virtual void SetVertex(const TVector3 &)
Definition: ParticleMC.cxx:299
Double32_t yv
y coordinate of particle production vertex
Definition: ParticleMC.h:385
virtual Double_t GetP() const
Definition: ParticleMC.h:460
Double32_t phi
Azimuthal angle.
Definition: ParticleMC.h:393
Double32_t pt
Transverse momentum of particle.
Definition: ParticleMC.h:383
Double32_t rapidity
Rapidity of particle.
Definition: ParticleMC.h:394
virtual const ParticleMC * GetParent() const
Definition: ParticleMC.cxx:205
virtual Double_t GetPhi() const
Definition: ParticleMC.h:468
virtual Double_t GetPy() const
Definition: ParticleMC.h:436
virtual ~ParticleMC()
Definition: ParticleMC.cxx:100
virtual Double_t GetXFeynman() const
Definition: ParticleMC.h:484
virtual void Set4Vector(const TLorentzVector &)
Definition: ParticleMC.cxx:284
Double32_t theta
Polar angle.
Definition: ParticleMC.h:392
virtual const ParticleMC * GetChild(UShort_t) const
Definition: ParticleMC.cxx:184
virtual Pid GetPdgCode() const
Definition: ParticleMC.h:272
virtual void SetStatus(int i)
Definition: ParticleMC.h:314
virtual Bool_t HasChild(Int_t) const
Definition: ParticleMC.cxx:217
virtual void ComputeEventDependentQuantities(EventMC &)
Definition: ParticleMC.cxx:132
virtual TLorentzVector Get4VectorInHadronBosonFrame() const
Definition: ParticleMC.cxx:229
Double32_t phiPrf
Definition: ParticleMC.h:403
UShort_t daughter
I of first child particle.
Definition: ParticleMC.h:355
Int_t id
PDG particle code.
Definition: ParticleMC.h:353
virtual Double_t GetEta() const
Definition: ParticleMC.h:476
virtual UShort_t GetChildNIndex() const
Definition: ParticleMC.h:428
virtual UInt_t GetNChildren() const
Definition: ParticleMC.h:500
Double32_t p
Total momentum of particle.
Definition: ParticleMC.h:391
virtual Pid Id() const
Definition: ParticleMC.h:496
Int_t parentId
PDG code of this particle's parent.
Definition: ParticleMC.h:389
virtual UInt_t GetIndex() const
Definition: ParticleMC.h:412
virtual Double_t GetPt() const
Definition: ParticleMC.h:448
virtual Double_t GetTheta() const
Definition: ParticleMC.h:464
virtual Double_t GetE() const
Definition: ParticleMC.h:506
Definition: Pid.h:22
UShort_t KS
Particle status code: see PYTHIA manual.
Definition: ParticleMC.h:352
Double32_t xFeynman
Feynman x = pz/(2sqrt(s))
Definition: ParticleMC.h:398
virtual void SetParentIndex(int i)
Definition: ParticleMC.h:323
Abstract base class for a general particle.
virtual void SetChild1Index(int i)
Definition: ParticleMC.h:327
virtual void SetIndex(int i)
Definition: ParticleMC.h:310
virtual Double_t GetPx() const
Definition: ParticleMC.h:432
Double32_t z
Definition: ParticleMC.h:396
virtual void SetParentId(int i)
Definition: ParticleMC.h:347
virtual TLorentzVector PxPyPzE() const
Definition: ParticleMC.h:222
virtual UShort_t GetParentIndex() const
Definition: ParticleMC.h:420