Belle II Software  release-05-01-25
SensitiveDetectorBase.h
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2010 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Martin Ritter *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #ifndef VXD_SENSITIVEDETECTORBASE_H
12 #define VXD_SENSITIVEDETECTORBASE_H
13 
14 //If this is defined, we create a root file with all information necessary to
15 //check functionality of the implementation
16 //#define VXD_SENSITIVEDETECTOR_DEBUG
17 
18 #include <simulation/kernel/SensitiveDetectorBase.h>
19 #include <vxd/simulation/SensorTraversal.h>
20 #include <vxd/geometry/SensorInfoBase.h>
21 #include <stack>
22 #include <vector>
23 
24 namespace Belle2 {
30  namespace VXD {
39  class SensitiveDetectorBase: public Simulation::SensitiveDetectorBase {
40  public:
48  Simulation::SensitiveDetectorBase((info->getType() == SensorInfoBase::PXD ? "PXD " : "SVD ") + (std::string)info->getID(),
49  info->getType() == SensorInfoBase::PXD ? Const::PXD : Const::SVD), m_info(info) {};
50 
62  void setOptions(bool seeNeutrons, bool onlyPrimaryTrueHits, float distanceTolerance,
63  float electronTolerance, float minimumElectrons)
64  {
65  m_seeNeutrons = seeNeutrons;
66  m_onlyPrimaryTrueHits = onlyPrimaryTrueHits;
67  m_distanceTolerance = distanceTolerance;
68  m_electronTolerance = electronTolerance;
69  m_minimumElectrons = minimumElectrons;
70  }
71 
73  virtual ~SensitiveDetectorBase() override
74  {
75  if (m_info) delete m_info;
76  }
77 
79  SensorInfoBase* getSensorInfo() { return m_info; }
80 
82  VxdID getSensorID() const { return m_info->getID(); }
83 
84  protected:
92  virtual int saveTrueHit(const SensorTraversal& traversal) = 0;
93 
103  virtual int saveSimHit(const SensorTraversal& traversal, const SensorTraversal::range& points) = 0;
104 
112  virtual void saveRelations(const SensorTraversal& traversal, int trueHitIndex,
113  std::vector<std::pair<unsigned int, float>> simHitIndices) = 0;
114 
127  std::vector<unsigned int> simplifyEnergyDeposit(const SensorTraversal::range& points);
128 
139  StepInformation findMidPoint(const SensorTraversal& traversal);
140  private:
146  bool step(G4Step* step, G4TouchableHistory*) override;
147 
151  bool finishTrack();
152 
162  std::vector<std::pair<unsigned int, float>> createSimHits();
163 
165  SensorInfoBase* m_info {0};
167  std::stack<SensorTraversal> m_tracks;
170  float m_distanceTolerance {0};
176  float m_minimumElectrons {0};
179  bool m_seeNeutrons {false};
181  bool m_onlyPrimaryTrueHits {false};
182  };
183  }
185 } //Belle2 namespace
186 #endif
Belle2::VXD::SensitiveDetectorBase::m_minimumElectrons
float m_minimumElectrons
minimum number of electrons a track must deposit for SimHit/TrueHits to be created
Definition: SensitiveDetectorBase.h:184
Belle2::VXD::SensitiveDetectorBase::m_info
SensorInfoBase * m_info
Pointer to the SensorInfo associated with this instance.
Definition: SensitiveDetectorBase.h:173
Belle2::VXD::SensitiveDetectorBase::~SensitiveDetectorBase
virtual ~SensitiveDetectorBase() override
Destructor freeing the sensor Info.
Definition: SensitiveDetectorBase.h:81
Belle2::getID
int getID(const std::vector< double > &breaks, double t)
get id of the time point t
Definition: calibTools.h:71
Belle2::VXD::SensitiveDetectorBase::createSimHits
std::vector< std::pair< unsigned int, float > > createSimHits()
Determine which SimHits to create.
Definition: SensitiveDetectorBase.cc:133
Belle2::VXD::SensitiveDetectorBase::finishTrack
bool finishTrack()
Process a track once all steps are known.
Definition: SensitiveDetectorBase.cc:102
Belle2::VXD::SensitiveDetectorBase::getSensorID
VxdID getSensorID() const
Return the VxdID belonging to this sensitive detector.
Definition: SensitiveDetectorBase.h:90
Belle2::VXD::SensitiveDetectorBase
Base class for Sensitive Detector implementation of PXD and SVD.
Definition: SensitiveDetectorBase.h:47
Belle2::VXD::SensorInfoBase
Base class to provide Sensor Information for PXD and SVD.
Definition: SensorInfoBase.h:40
Belle2::VXD::SensorInfoBase::getID
VxdID getID() const
Return the ID of the Sensor.
Definition: SensorInfoBase.h:85
Belle2::VXD::SensitiveDetectorBase::SensitiveDetectorBase
SensitiveDetectorBase(SensorInfoBase *info)
Constructor.
Definition: SensitiveDetectorBase.h:55
Belle2::VXD::SensitiveDetectorBase::findMidPoint
StepInformation findMidPoint(const SensorTraversal &traversal)
Find the mid-point of the track traversal.
Definition: SensitiveDetectorBase.cc:281
Belle2::VXD::SensitiveDetectorBase::m_seeNeutrons
bool m_seeNeutrons
also create SimHit/TrueHit objects for neutrons (or charged particles which deposit less than m_minim...
Definition: SensitiveDetectorBase.h:187
Belle2::VXD::SensitiveDetectorBase::m_electronTolerance
float m_electronTolerance
maximum relative difference between electron density of two steps where they can be considered simila...
Definition: SensitiveDetectorBase.h:181
Belle2::VXD::SensitiveDetectorBase::m_onlyPrimaryTrueHits
bool m_onlyPrimaryTrueHits
only create TrueHits for primary particles
Definition: SensitiveDetectorBase.h:189
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::SensorTraversal::range
std::pair< iterator, iterator > range
Iterator pair for a set of points.
Definition: SensorTraversal.h:57
Belle2::VXD::SensitiveDetectorBase::step
bool step(G4Step *step, G4TouchableHistory *) override
Process a single Geant4 Step.
Definition: SensitiveDetectorBase.cc:37
Belle2::VXD::SensitiveDetectorBase::setOptions
void setOptions(bool seeNeutrons, bool onlyPrimaryTrueHits, float distanceTolerance, float electronTolerance, float minimumElectrons)
Set all common options.
Definition: SensitiveDetectorBase.h:70
Belle2::VXD::SensitiveDetectorBase::saveTrueHit
virtual int saveTrueHit(const SensorTraversal &traversal)=0
Save the actual TrueHit for this sensor traversal.
Belle2::VXD::SensitiveDetectorBase::getSensorInfo
SensorInfoBase * getSensorInfo()
Return a pointer to the SensorInfo associated with this instance.
Definition: SensitiveDetectorBase.h:87
Belle2::VXD::SensitiveDetectorBase::simplifyEnergyDeposit
std::vector< unsigned int > simplifyEnergyDeposit(const SensorTraversal::range &points)
Simplify the energy deposition profile using Douglas-Peuker-Algorithm We normally force a Geant4 step...
Definition: SensitiveDetectorBase.cc:190
Belle2::VXD::SensitiveDetectorBase::m_tracks
std::stack< SensorTraversal > m_tracks
stack of SensorTraversal information for all tracks not finished so far
Definition: SensitiveDetectorBase.h:175
Belle2::VXD::SensitiveDetectorBase::saveSimHit
virtual int saveSimHit(const SensorTraversal &traversal, const SensorTraversal::range &points)=0
Save a SimHit for this track including the given points.
Belle2::Const
This class provides a set of constants for the framework.
Definition: Const.h:36
Belle2::VXD::SensitiveDetectorBase::saveRelations
virtual void saveRelations(const SensorTraversal &traversal, int trueHitIndex, std::vector< std::pair< unsigned int, float >> simHitIndices)=0
Save the relations between MCParticle, TrueHit and SimHits.
Belle2::VXD::SensitiveDetectorBase::m_distanceTolerance
float m_distanceTolerance
maximum distance between step point and linear interpolation of sensor traversal before a new simhit ...
Definition: SensitiveDetectorBase.h:178