Belle II Software development
SensitiveDetector.h
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9#pragma once
10
11/* KLM headers. */
12#include <klm/bklm/geometry/GeometryPar.h>
13#include <klm/bklm/geometry/Module.h>
14#include <klm/dataobjects/eklm/EKLMElementNumbers.h>
15#include <klm/dataobjects/KLMSimHit.h>
16#include <klm/dbobjects/bklm/BKLMSimulationPar.h>
17
18/* Basf2 headers. */
19#include <framework/datastore/StoreArray.h>
20#include <framework/datastore/RelationArray.h>
21#include <framework/database/DBObjPtr.h>
22#include <mdst/dataobjects/MCParticle.h>
23#include <simulation/kernel/SensitiveDetectorBase.h>
24
25namespace Belle2 {
31 class BkgSensitiveDetector;
32
33 namespace KLM {
34
39
40 public:
41
47 SensitiveDetector(const G4String& name, KLMSubdetectorNumber subdetector);
48
52 bool step(G4Step*, G4TouchableHistory*) override;
53
54 private:
55
59 bool stepBKLM(G4Step*, G4TouchableHistory*);
60
64 bool stepEKLM(G4Step*, G4TouchableHistory*);
65
69 void convertHitToRPCStrips(const CLHEP::Hep3Vector&, const bklm::Module*, int&, int&, int&, int&);
70
73
76
79
81 static constexpr int m_DepthSection = 2;
82
84 static constexpr int m_DepthSector = 3;
85
87 static constexpr int m_DepthLayer = 5;
88
90 static constexpr int m_DepthPlane = 9;
91
93 static constexpr int m_DepthScintillator = 10;
94
100
103
106
109
112
115
118
119 };
120
121 }
122
124}
The Class for BeamBackground Sensitive Detector.
Class for accessing objects in the database.
Definition: DBObjPtr.h:21
EKLM element numbers.
KLM sensitive-detector class.
bool stepBKLM(G4Step *, G4TouchableHistory *)
Step in BKLM.
DBObjPtr< BKLMSimulationPar > m_SimPar
Simulation parameters (from DB).
bool stepEKLM(G4Step *, G4TouchableHistory *)
Step in EKLM.
bool m_FirstCall
Flag to enforce once-only initializations in Initialize().
KLMSubdetectorNumber m_Subdetector
Subdetector.
RelationArray m_MCParticlesToKLMSimHits
Relation array between MCPartices and KLMSimHits.
static constexpr int m_DepthSection
Section depth.
static constexpr int m_DepthSector
Sector depth.
bklm::GeometryPar * m_GeoPar
Pointer to GeometryPar singleton.
const EKLMElementNumbers * m_ElementNumbers
EKLM element numbers.
bool step(G4Step *, G4TouchableHistory *) override
Process each step and store KLMSimHits.
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.
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.
StoreArray< MCParticle > m_MCParticles
MC particles.
static constexpr int m_DepthLayer
Layer depth.
StoreArray< KLMSimHit > m_KLMSimHits
Simulation hits.
Low-level class to create/modify relations between StoreArrays.
Definition: RelationArray.h:62
Base class for all Sensitive Detectors to create hits during simulation.
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
Provides BKLM geometry parameters for simulation, reconstruction etc (from Gearbox or DataBase)
Definition: GeometryPar.h:37
Define the geometry of a BKLM module Each sector [octant] contains Modules.
Definition: Module.h:76
uint16_t KLMSubdetectorNumber
Subdetector number.
Abstract base class for different kinds of events.