Belle II Software development
SensitiveDetector Class Reference

KLM sensitive-detector class. More...

#include <SensitiveDetector.h>

Inheritance diagram for SensitiveDetector:
SensitiveDetectorBase

Public Member Functions

 SensitiveDetector (const G4String &name, KLMSubdetectorNumber subdetector)
 Constructor.
 
bool step (G4Step *, G4TouchableHistory *) override
 Process each step and store KLMSimHits.
 

Static Public Member Functions

static const std::map< std::string, RelationArray::EConsolidationAction > & getMCParticleRelations ()
 Return a list of all registered Relations with MCParticles.
 
static void setActive (bool activeStatus)
 Enable/Disable all Sensitive Detectors.
 
static void registerMCParticleRelation (const std::string &name, RelationArray::EConsolidationAction ignoreAction=RelationArray::c_negativeWeight)
 Register an relation involving MCParticles.
 
static void registerMCParticleRelation (const RelationArray &relation, RelationArray::EConsolidationAction ignoreAction=RelationArray::c_negativeWeight)
 Overload to make it easer to register MCParticle relations.
 

Private Member Functions

bool stepBKLM (G4Step *, G4TouchableHistory *)
 Step in BKLM.
 
bool stepEKLM (G4Step *, G4TouchableHistory *)
 Step in EKLM.
 
void convertHitToRPCStrips (const CLHEP::Hep3Vector &, const bklm::Module *, int &, int &, int &, int &)
 Find the ranges of matching RPC strips for BKLM hits.
 
virtual bool ProcessHits (G4Step *aStep, G4TouchableHistory *aROhist)
 Check if recording hits is enabled and if so call step() and set the correct MCParticle flag.
 

Private Attributes

KLMSubdetectorNumber m_Subdetector
 Subdetector.
 
bool m_FirstCall
 Flag to enforce once-only initializations in Initialize().
 
double m_HitTimeMax
 Maximum permissible hit time based on overflow of LeCroy 1877 TDC.
 
BkgSensitiveDetectorm_BkgSensitiveDetector
 Pointer to a sensitive-detector object used for beam-background steps.
 
bklm::GeometryParm_GeoPar
 Pointer to GeometryPar singleton.
 
const EKLMElementNumbersm_ElementNumbers
 EKLM element numbers.
 
DBObjPtr< BKLMSimulationParm_SimPar
 Simulation parameters (from DB).
 
StoreArray< MCParticlem_MCParticles
 MC particles.
 
StoreArray< KLMSimHitm_KLMSimHits
 Simulation hits.
 
RelationArray m_MCParticlesToKLMSimHits {m_MCParticles, m_KLMSimHits}
 Relation array between MCPartices and KLMSimHits.
 
Const::EDetector m_subdetector
 Subdetector the class belongs to.
 

Static Private Attributes

static constexpr int m_DepthSection = 2
 Section depth.
 
static constexpr int m_DepthSector = 3
 Sector depth.
 
static constexpr int m_DepthLayer = 5
 Layer depth.
 
static constexpr int m_DepthPlane = 9
 Plane depth.
 
static constexpr int m_DepthScintillator = 10
 Scintillator depth.
 
static std::map< std::string, RelationArray::EConsolidationActions_mcRelations
 Static set holding all relations which have to be updated at the end of the Event.
 
static bool s_active
 Static bool which indicates wether recording of hits is enabled.
 

Detailed Description

KLM sensitive-detector class.

Definition at line 38 of file SensitiveDetector.h.

Constructor & Destructor Documentation

◆ SensitiveDetector()

SensitiveDetector ( const G4String &  name,
KLMSubdetectorNumber  subdetector 
)

Constructor.

Parameters
[in]nameSensitive-detector name.
[in]subdetectorSubdetector.

Definition at line 33 of file SensitiveDetector.cc.

34 :
35 SensitiveDetectorBase(name, Const::KLM),
36 m_Subdetector(subdetector),
37 m_FirstCall(true),
39 m_GeoPar(nullptr)
40{
41 if (!m_SimPar.isValid())
42 B2FATAL("BKLM simulation parameters are not available.");
43 m_HitTimeMax = m_SimPar->getHitTimeMax();
45 m_KLMSimHits.registerInDataStore();
48}
DBObjPtr< BKLMSimulationPar > m_SimPar
Simulation parameters (from DB).
bool m_FirstCall
Flag to enforce once-only initializations in Initialize().
KLMSubdetectorNumber m_Subdetector
Subdetector.
RelationArray m_MCParticlesToKLMSimHits
Relation array between MCPartices and KLMSimHits.
bklm::GeometryPar * m_GeoPar
Pointer to GeometryPar singleton.
double m_HitTimeMax
Maximum permissible hit time based on overflow of LeCroy 1877 TDC.
BkgSensitiveDetector * m_BkgSensitiveDetector
Pointer to a sensitive-detector object used for beam-background steps.
StoreArray< MCParticle > m_MCParticles
MC particles.
StoreArray< KLMSimHit > m_KLMSimHits
Simulation hits.
static void registerMCParticleRelation(const std::string &name, RelationArray::EConsolidationAction ignoreAction=RelationArray::c_negativeWeight)
Register an relation involving MCParticles.
SensitiveDetectorBase(const std::string &name, Const::EDetector subdetector)
Create a new Sensitive detecor with a given name and belonging to a given subdetector.
bool isOptional(const std::string &name="")
Tell the DataStore about an optional input.
bool registerInDataStore(DataStore::EStoreFlags storeFlags=DataStore::c_WriteOut)
Register the object/array in the DataStore.

Member Function Documentation

◆ convertHitToRPCStrips()

void convertHitToRPCStrips ( const CLHEP::Hep3Vector &  localPosition,
const bklm::Module m,
int &  phiStripLower,
int &  phiStripUpper,
int &  zStripLower,
int &  zStripUpper 
)
private

Find the ranges of matching RPC strips for BKLM hits.

Definition at line 239 of file SensitiveDetector.cc.

242{
243 double phiStripD = m->getPhiStrip(localPosition);
244 int phiStrip = int(phiStripD);
245 int pMin = m->getPhiStripMin();
246 if (phiStrip < pMin)
247 return;
248 int pMax = m->getPhiStripMax();
249 if (phiStrip > pMax)
250 return;
251
252 double zStripD = m->getZStrip(localPosition);
253 int zStrip = int(zStripD);
254 int zMin = m->getZStripMin();
255 if (zStrip < zMin)
256 return;
257 int zMax = m->getZStripMax();
258 if (zStrip > zMax)
259 return;
260
261 phiStripLower = phiStrip;
262 phiStripUpper = phiStrip;
263 zStripLower = zStrip;
264 zStripUpper = zStrip;
265 double phiStripDiv = fmod(phiStripD, 1.0) - 0.5; // between -0.5 and +0.5 within central phiStrip
266 double zStripDiv = fmod(zStripD, 1.0) - 0.5; // between -0.5 and +0.5 within central zStrip
267 int n = 0;
268 double rand = gRandom->Uniform();
269 for (n = 1; n < m_SimPar->getMaxMultiplicity(); ++n) {
270 if (m_SimPar->getPhiMultiplicityCDF(phiStripDiv, n) > rand)
271 break;
272 }
273 int nextStrip = (phiStripDiv > 0.0 ? 1 : -1);
274 while (--n > 0) {
275 phiStrip += nextStrip;
276 if ((phiStrip >= pMin) && (phiStrip <= pMax)) {
277 phiStripLower = std::min(phiStrip, phiStripLower);
278 phiStripUpper = std::max(phiStrip, phiStripUpper);
279 }
280 nextStrip = (nextStrip > 0 ? -(1 + nextStrip) : 1 - nextStrip);
281 }
282 rand = gRandom->Uniform();
283 for (n = 1; n < m_SimPar->getMaxMultiplicity(); ++n) {
284 if (m_SimPar->getZMultiplicityCDF(zStripDiv, n) > rand)
285 break;
286 }
287 nextStrip = (zStripDiv > 0.0 ? 1 : -1);
288 while (--n > 0) {
289 zStrip += nextStrip;
290 if ((zStrip >= zMin) && (zStrip <= zMax)) {
291 zStripLower = std::min(zStrip, zStripLower);
292 zStripUpper = std::max(zStrip, zStripUpper);
293 }
294 nextStrip = (nextStrip > 0 ? -(1 + nextStrip) : 1 - nextStrip);
295 }
296 return;
297}

◆ getMCParticleRelations()

static const std::map< std::string, RelationArray::EConsolidationAction > & getMCParticleRelations ( )
inlinestaticinherited

Return a list of all registered Relations with MCParticles.

Definition at line 42 of file SensitiveDetectorBase.h.

42{ return s_mcRelations; }
static std::map< std::string, RelationArray::EConsolidationAction > s_mcRelations
Static set holding all relations which have to be updated at the end of the Event.

◆ ProcessHits()

bool ProcessHits ( G4Step *  aStep,
G4TouchableHistory *  aROhist 
)
inlineprivatevirtualinherited

Check if recording hits is enabled and if so call step() and set the correct MCParticle flag.

Called by Geant4 for each step inside the sensitive volumes attached

Definition at line 94 of file SensitiveDetectorBase.h.

95 {
96 if (!s_active) return false;
97 bool result = step(aStep, aROhist);
98 // Do not include hits from invalid detector (beast,teastbeam, etc.)
99 if (result && (m_subdetector != Const::invalidDetector)) TrackInfo::getInfo(*aStep->GetTrack()).addSeenInDetector(m_subdetector);
100 return result;
101 }
virtual bool step(G4Step *step, G4TouchableHistory *ROhist)=0
Process a Geant4 step in any of the sensitive volumes attached to this sensitive detector.
Const::EDetector m_subdetector
Subdetector the class belongs to.
static bool s_active
Static bool which indicates wether recording of hits is enabled.
static Payload getInfo(Carrier &obj)
Static function to just return UserInformation attached to the obj of type Carrier.
Definition: UserInfo.h:100

◆ registerMCParticleRelation() [1/2]

static void registerMCParticleRelation ( const RelationArray relation,
RelationArray::EConsolidationAction  ignoreAction = RelationArray::c_negativeWeight 
)
inlinestaticinherited

Overload to make it easer to register MCParticle relations.

Parameters
relationRelationArray to register
ignoreAction

Definition at line 66 of file SensitiveDetectorBase.h.

67 { registerMCParticleRelation(relation.getName(), ignoreAction); }

◆ registerMCParticleRelation() [2/2]

void registerMCParticleRelation ( const std::string &  name,
RelationArray::EConsolidationAction  ignoreAction = RelationArray::c_negativeWeight 
)
staticinherited

Register an relation involving MCParticles.

All Relations which point from an MCParticle to something have to be registered with addMCParticleRelation() because the index of the MCParticles might change at the end of the event. During simulation, the TrackID should be used as index of the MCParticle

Parameters
nameName of the relation to register
ignoreAction

Definition at line 22 of file SensitiveDetectorBase.cc.

23 {
24 std::pair<std::map<std::string, RelationArray::EConsolidationAction>::iterator, bool> insert = s_mcRelations.insert(std::make_pair(
25 name, ignoreAction));
26 //If the relation already exists and the ignoreAction is different we do have a problem
27 if (!insert.second && insert.first->second != ignoreAction) {
28 B2FATAL("MCParticle Relation " << name << " already registered with different ignore action.");
29 }
30 }

◆ setActive()

static void setActive ( bool  activeStatus)
inlinestaticinherited

Enable/Disable all Sensitive Detectors.

By default, all sensitive detectors won't create hits to make it possible to use the Geant4 Navigator for non-simulation purposes. Only during simulation the sensitive detectors will be enabled to record hits

Parameters
activeStatusbool to indicate wether hits should be recorded

Definition at line 50 of file SensitiveDetectorBase.h.

50{ s_active = activeStatus; }

◆ step()

bool step ( G4Step *  step,
G4TouchableHistory *  history 
)
overridevirtual

Process each step and store KLMSimHits.

Implements SensitiveDetectorBase.

Definition at line 299 of file SensitiveDetector.cc.

300{
302 return stepBKLM(step, history);
303 else
304 return stepEKLM(step, history);
305}
bool stepBKLM(G4Step *, G4TouchableHistory *)
Step in BKLM.
bool stepEKLM(G4Step *, G4TouchableHistory *)
Step in EKLM.
bool step(G4Step *, G4TouchableHistory *) override
Process each step and store KLMSimHits.

◆ stepBKLM()

G4bool stepBKLM ( G4Step *  step,
G4TouchableHistory *  history 
)
private

Step in BKLM.

Definition at line 117 of file SensitiveDetector.cc.

118{
119 /* Once-only initializations (constructor is called too early for these). */
120 if (m_FirstCall) {
121 m_FirstCall = false;
124 m_BkgSensitiveDetector = new BkgSensitiveDetector("BKLM");
125 if (!gRandom)
126 B2FATAL("gRandom is not initialized; please set up gRandom first");
127 }
128
129 // Record a BeamBackHit for any particle
130 if (m_BkgSensitiveDetector != nullptr) {
132 }
133
134 // It is not necessary to detect motion from one volume to another (or track death
135 // in the RPC gas volume). Experimentation shows that most tracks pass through the
136 // RPC gas volume in one step (although, on occasion, a delta ray will take a couple
137 // of short steps entirely within gas). Therefore, save every step in the gas
138 // instead of trying to find entry and exit points and then saving only the midpoint.
139 // Do same for scintillators.
140
141 double eDep = step->GetTotalEnergyDeposit() / CLHEP::MeV; // GEANT4: in MeV
142 G4StepPoint* preStep = step->GetPreStepPoint();
143 G4StepPoint* postStep = step->GetPostStepPoint();
144 G4Track* track = step->GetTrack();
145
146 // Record a KLMSimHit for a charged track that deposits some energy.
147 if ((eDep > 0.0) && (postStep->GetCharge() != 0.0)) {
148 const G4VTouchable* hist = preStep->GetTouchable();
149 int depth = hist->GetHistoryDepth();
150 if (depth < m_DepthPlane) {
151 B2WARNING("BKLM SensitiveDetector::step(): "
152 << LogVar("Touchable HistoryDepth", depth)
153 << LogVar("Should be at least", m_DepthPlane));
154 return false;
155 }
156 int plane = hist->GetCopyNumber(depth - m_DepthPlane);
157 int layer = hist->GetCopyNumber(depth - m_DepthLayer);
158 int sector = hist->GetCopyNumber(depth - m_DepthSector);
159 int section = hist->GetCopyNumber(depth - m_DepthSection);
160 double time = 0.5 * (preStep->GetGlobalTime() + postStep->GetGlobalTime()); // GEANT4: in ns
161 if (time > m_HitTimeMax)
162 return false;
163 const CLHEP::Hep3Vector globalPosition = 0.5 * (preStep->GetPosition() + postStep->GetPosition()) / CLHEP::cm; // in cm
164 const bklm::Module* m = m_GeoPar->findModule(section, sector, layer);
165 const CLHEP::Hep3Vector localPosition = m->globalToLocal(globalPosition);
166 int trackID = track->GetTrackID();
167 if (m->hasRPCs()) {
168 const CLHEP::Hep3Vector propagationTimes =
169 m->getPropagationTimes(localPosition);
170 int phiStripLower = -1;
171 int phiStripUpper = -1;
172 int zStripLower = -1;
173 int zStripUpper = -1;
174 convertHitToRPCStrips(localPosition, m, phiStripLower, phiStripUpper, zStripLower, zStripUpper);
175 if (zStripLower > 0) {
176 KLMSimHit* simHit = m_KLMSimHits.appendNew();
177 simHit->setSubdetector(KLMElementNumbers::c_BKLM);
178 simHit->setSection(section);
179 simHit->setSector(sector);
180 simHit->setLayer(layer);
181 simHit->setPlane(BKLMElementNumbers::c_ZPlane);
182 simHit->setStrip(zStripLower);
183 simHit->setLastStrip(zStripUpper);
184 simHit->setPropagationTime(propagationTimes.z());
185 simHit->setTime(time);
186 simHit->setEnergyDeposit(eDep);
187 simHit->setPosition(globalPosition.x(), globalPosition.y(),
188 globalPosition.z());
189 simHit->setPDG(track->GetDefinition()->GetPDGEncoding());
190 m_MCParticlesToKLMSimHits.add(trackID, simHit->getArrayIndex());
191 }
192 if (phiStripLower > 0) {
193 KLMSimHit* simHit = m_KLMSimHits.appendNew();
194 simHit->setSubdetector(KLMElementNumbers::c_BKLM);
195 simHit->setSection(section);
196 simHit->setSector(sector);
197 simHit->setLayer(layer);
198 simHit->setPlane(BKLMElementNumbers::c_PhiPlane);
199 simHit->setStrip(phiStripLower);
200 simHit->setLastStrip(phiStripUpper);
201 simHit->setPropagationTime(propagationTimes.y());
202 simHit->setTime(time);
203 simHit->setEnergyDeposit(eDep);
204 simHit->setPosition(globalPosition.x(), globalPosition.y(),
205 globalPosition.z());
206 simHit->setPDG(track->GetDefinition()->GetPDGEncoding());
207 m_MCParticlesToKLMSimHits.add(trackID, simHit->getArrayIndex());
208 }
209 } else {
210 int scint = hist->GetCopyNumber(depth - m_DepthScintillator);
211 bool phiPlane = (plane == BKLM_INNER);
212 double propagationTime =
213 m->getPropagationTime(localPosition, scint, phiPlane);
214 KLMSimHit* simHit = m_KLMSimHits.appendNew();
215 simHit->setSubdetector(KLMElementNumbers::c_BKLM);
216 simHit->setSection(section);
217 simHit->setSector(sector);
218 simHit->setLayer(layer);
219 if (phiPlane) {
220 simHit->setPlane(BKLMElementNumbers::c_PhiPlane);
221 } else {
222 simHit->setPlane(BKLMElementNumbers::c_ZPlane);
223 }
224 simHit->setStrip(scint);
225 simHit->setLastStrip(scint);
226 simHit->setPropagationTime(propagationTime);
227 simHit->setTime(time);
228 simHit->setEnergyDeposit(eDep);
229 simHit->setPosition(globalPosition.x(), globalPosition.y(),
230 globalPosition.z());
231 simHit->setPDG(track->GetDefinition()->GetPDGEncoding());
232 m_MCParticlesToKLMSimHits.add(trackID, simHit->getArrayIndex());
233 }
234 return true;
235 }
236 return false;
237}
static constexpr int m_DepthSection
Section depth.
static constexpr int m_DepthSector
Sector depth.
static constexpr int m_DepthScintillator
Scintillator depth.
static constexpr int m_DepthPlane
Plane depth.
void convertHitToRPCStrips(const CLHEP::Hep3Vector &, const bklm::Module *, int &, int &, int &, int &)
Find the ranges of matching RPC strips for BKLM hits.
static constexpr int m_DepthLayer
Layer depth.
void add(index_type from, index_type to, weight_type weight=1.0)
Add a new element to the relation.
bool doBeamBackgroundStudy(void) const
Get the beam background study flag.
Definition: GeometryPar.h:51
const Module * findModule(int section, int sector, int layer) const
Get the pointer to the definition of a module.
Definition: GeometryPar.cc:721
static GeometryPar * instance(void)
Static method to get a reference to the singleton GeometryPar instance.
Definition: GeometryPar.cc:27
Class to store variables with their name which were sent to the logging service.
bool step(G4Step *aStep, G4TouchableHistory *) override
Process each step and calculate variables defined in PXDSimHit.

◆ stepEKLM()

bool stepEKLM ( G4Step *  aStep,
G4TouchableHistory *  history 
)
private

Step in EKLM.

Definition at line 50 of file SensitiveDetector.cc.

51{
52 /* Once-only initializations (constructor is called too early for these). */
53 if (m_FirstCall) {
54 m_FirstCall = false;
55 const EKLM::GeometryData* geometryData = &EKLM::GeometryData::Instance();
56 if (geometryData->beamBackgroundStudy())
57 m_BkgSensitiveDetector = new BkgSensitiveDetector("EKLM");
58 }
59
60 // Record a BeamBackHit for any particle
61 if (m_BkgSensitiveDetector != nullptr) {
62 m_BkgSensitiveDetector->step(aStep, history);
63 }
64
65 const int stripLevel = 1;
66 int section, layer, sector, plane, strip, stripGlobal;
67 HepGeom::Point3D<double> gpos, lpos;
68 G4TouchableHandle hist = aStep->GetPreStepPoint()->GetTouchableHandle();
69 section = hist->GetVolume(stripLevel + 6)->GetCopyNo();
70 layer = hist->GetVolume(stripLevel + 5)->GetCopyNo();
71 sector = hist->GetVolume(stripLevel + 4)->GetCopyNo();
72 plane = hist->GetVolume(stripLevel + 3)->GetCopyNo();
73 strip = hist->GetVolume(stripLevel)->GetCopyNo();
74 stripGlobal = m_ElementNumbers->stripNumber(
75 section, layer, sector, plane, strip);
76 const G4double eDep = aStep->GetTotalEnergyDeposit();
77 /* Do not record hits without deposited energy. */
78 if (eDep <= 0)
79 return false;
80 const G4Track& track = * aStep->GetTrack();
81 const G4double hitTime = track.GetGlobalTime();
82 if (hitTime > m_HitTimeMax)
83 return false;
84 /* Hit position. */
85 gpos = 0.5 * (aStep->GetPostStepPoint()->GetPosition() +
86 aStep->GetPreStepPoint()->GetPosition());
87 lpos = hist->GetHistory()->GetTopTransform().TransformPoint(gpos);
88 /* Create step hit and store in to DataStore */
89 KLMSimHit* hit = m_KLMSimHits.appendNew();
90 hit->setSubdetector(KLMElementNumbers::c_EKLM);
91 hit->setSection(section);
92 hit->setSector(sector);
93 hit->setLayer(layer);
94 hit->setPlane(plane);
95 hit->setStrip(strip);
96 hit->setLastStrip(strip);
97 // hit->setPropagationTime();
98 hit->setTime(hitTime);
99 hit->setEnergyDeposit(eDep);
100 CLHEP::Hep3Vector trackMomentum = track.GetMomentum();
101 hit->setMomentum(ROOT::Math::PxPyPzEVector(trackMomentum.x(), trackMomentum.y(),
102 trackMomentum.z(), track.GetTotalEnergy()));
103 hit->setTrackID(track.GetTrackID());
104 hit->setParentTrackID(track.GetParentID());
105 hit->setLocalPosition(lpos.x() / CLHEP::mm * Unit::mm,
106 lpos.y() / CLHEP::mm * Unit::mm,
107 lpos.z() / CLHEP::mm * Unit::mm);
108 hit->setPosition(gpos.x() / CLHEP::mm * Unit::mm,
109 gpos.y() / CLHEP::mm * Unit::mm,
110 gpos.z() / CLHEP::mm * Unit::mm);
111 hit->setPDG(track.GetDefinition()->GetPDGEncoding());
112 hit->setVolumeID(stripGlobal);
113 m_MCParticlesToKLMSimHits.add(track.GetTrackID(), hit->getArrayIndex());
114 return true;
115}
int stripNumber(int section, int layer, int sector, int plane, int strip) const
Get strip number.
static const GeometryData & Instance(enum DataSource dataSource=c_Database, const GearDir *gearDir=nullptr)
Instantiation.
Definition: GeometryData.cc:33
const EKLMElementNumbers * m_ElementNumbers
EKLM element numbers.
static const double mm
[millimeters]
Definition: Unit.h:70

Member Data Documentation

◆ m_BkgSensitiveDetector

BkgSensitiveDetector* m_BkgSensitiveDetector
private

Pointer to a sensitive-detector object used for beam-background steps.

Definition at line 99 of file SensitiveDetector.h.

◆ m_DepthLayer

constexpr int m_DepthLayer = 5
staticconstexprprivate

Layer depth.

Definition at line 87 of file SensitiveDetector.h.

◆ m_DepthPlane

constexpr int m_DepthPlane = 9
staticconstexprprivate

Plane depth.

Definition at line 90 of file SensitiveDetector.h.

◆ m_DepthScintillator

constexpr int m_DepthScintillator = 10
staticconstexprprivate

Scintillator depth.

Definition at line 93 of file SensitiveDetector.h.

◆ m_DepthSection

constexpr int m_DepthSection = 2
staticconstexprprivate

Section depth.

Definition at line 81 of file SensitiveDetector.h.

◆ m_DepthSector

constexpr int m_DepthSector = 3
staticconstexprprivate

Sector depth.

Definition at line 84 of file SensitiveDetector.h.

◆ m_ElementNumbers

const EKLMElementNumbers* m_ElementNumbers
private

EKLM element numbers.

Definition at line 105 of file SensitiveDetector.h.

◆ m_FirstCall

bool m_FirstCall
private

Flag to enforce once-only initializations in Initialize().

Definition at line 75 of file SensitiveDetector.h.

◆ m_GeoPar

bklm::GeometryPar* m_GeoPar
private

Pointer to GeometryPar singleton.

Definition at line 102 of file SensitiveDetector.h.

◆ m_HitTimeMax

double m_HitTimeMax
private

Maximum permissible hit time based on overflow of LeCroy 1877 TDC.

Definition at line 78 of file SensitiveDetector.h.

◆ m_KLMSimHits

StoreArray<KLMSimHit> m_KLMSimHits
private

Simulation hits.

Definition at line 114 of file SensitiveDetector.h.

◆ m_MCParticles

StoreArray<MCParticle> m_MCParticles
private

MC particles.

Definition at line 111 of file SensitiveDetector.h.

◆ m_MCParticlesToKLMSimHits

RelationArray m_MCParticlesToKLMSimHits {m_MCParticles, m_KLMSimHits}
private

Relation array between MCPartices and KLMSimHits.

Definition at line 117 of file SensitiveDetector.h.

◆ m_SimPar

DBObjPtr<BKLMSimulationPar> m_SimPar
private

Simulation parameters (from DB).

Definition at line 108 of file SensitiveDetector.h.

◆ m_Subdetector

KLMSubdetectorNumber m_Subdetector
private

Subdetector.

Definition at line 72 of file SensitiveDetector.h.

◆ m_subdetector

Const::EDetector m_subdetector
privateinherited

Subdetector the class belongs to.

Definition at line 91 of file SensitiveDetectorBase.h.

◆ s_active

bool s_active
staticprivateinherited

Static bool which indicates wether recording of hits is enabled.

Definition at line 89 of file SensitiveDetectorBase.h.

◆ s_mcRelations

map< string, RelationArray::EConsolidationAction > s_mcRelations
staticprivateinherited

Static set holding all relations which have to be updated at the end of the Event.

Definition at line 87 of file SensitiveDetectorBase.h.


The documentation for this class was generated from the following files: