ATLAS Offline Software
Loading...
Searching...
No Matches
LUCID_SimHit.h
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2021 CERN for the benefit of the ATLAS collaboration
3*/
4
5#ifndef LUCID_SIMHIT_H
6#define LUCID_SIMHIT_H 1
7
10#include "CLHEP/Vector/ThreeVector.h"
11#include <string>
12#include "CLHEP/Random/RandFlat.h"
13
14// LUCID_SimHit hold all information needed by the digitization algorithms
15// to construct raw data objects.
16
18
19 public:
20
23 LUCID_SimHit(short tubeID,
24 int pdgCode,
25 int truthBarcode,
26 int genVolume,
27 float stepStartPosX,
28 float stepStartPosY,
29 float stepStartPosZ,
30 float stepEndPosX,
31 float stepEndPosY,
32 float stepEndPosZ,
33 float preStepTime,
34 float postStepTime,
35 float wavelength,
36 float energy);
37 LUCID_SimHit(short tubeID,
38 int pdgCode,
39 const HepMcParticleLink& partLink,
40 int genVolume,
41 float stepStartPosX,
42 float stepStartPosY,
43 float stepStartPosZ,
44 float stepEndPosX,
45 float stepEndPosY,
46 float stepEndPosZ,
47 float preStepTime,
48 float postStepTime,
49 float wavelength,
50 float energy);
51
52 std::string print() const;
53
54 const LUCID_SimHit& operator=(const LUCID_SimHit& t);
55
56 int operator == (const LUCID_SimHit&) const;
57 bool operator < (const LUCID_SimHit&) const;
58
59 const HepMcParticleLink& particleLink() const;
60
61 inline float GetX() const { return m_stepStartPosX; }
62 inline float GetY() const { return m_stepStartPosY; }
63 inline float GetZ() const { return m_stepStartPosZ; }
64
65 inline float GetEPX() const { return m_stepEndPosX; }
66 inline float GetEPY() const { return m_stepEndPosY; }
67 inline float GetEPZ() const { return m_stepEndPosZ; }
68
69 inline short GetTubeID () const { return m_tubeID; }
70 inline int truthBarcode () const { return m_partLink.barcode(); }
71 inline int truthID () const { return m_partLink.id(); }
72 inline int GetPdgCode () const { return m_pdgCode; }
73 inline int GetGenVolume () const { return m_genVolume; }
74 inline float GetPreStepTime () const { return m_preStepTime; }
75 inline float GetPostStepTime() const { return m_postStepTime; }
76 inline float GetWavelength () const { return m_wavelength; }
77 inline float GetEnergy () const { return m_energy; }
78
79 bool isDetected(CLHEP::HepRandomEngine* rndEngine) const;
80
81 private:
82
83 short m_tubeID{};
84 int m_pdgCode{};
85 HepMcParticleLink m_partLink; // link to the particle generating the hit
95 float m_wavelength{};
96 float m_energy{};
97};
98
100
101#endif
float GetEPY() const
float m_preStepTime
bool isDetected(CLHEP::HepRandomEngine *rndEngine) const
float GetEnergy() const
float GetPreStepTime() const
const HepMcParticleLink & particleLink() const
const LUCID_SimHit & operator=(const LUCID_SimHit &t)
float GetEPZ() const
float GetX() const
int truthID() const
float m_stepStartPosX
int operator==(const LUCID_SimHit &) const
bool operator<(const LUCID_SimHit &) const
HepMcParticleLink m_partLink
float m_postStepTime
float GetPostStepTime() const
float m_stepStartPosZ
float GetY() const
float m_stepEndPosZ
float GetEPX() const
int truthBarcode() const
int GetPdgCode() const
float m_stepEndPosX
float m_wavelength
float GetWavelength() const
short GetTubeID() const
float m_stepStartPosY
std::string print() const
float GetZ() const
float m_stepEndPosY
int GetGenVolume() const