Belle II Software development
SensitiveDetector Class Reference

Sensitive Detector implementation of the FANGS detector. More...

#include <SensitiveDetector.h>

Inheritance diagram for SensitiveDetector:
SensitiveDetectorBase

Public Member Functions

 SensitiveDetector ()
 Constructor.
 

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.
 

Protected Member Functions

bool step (G4Step *step, G4TouchableHistory *) override
 Step processing method.
 

Private Member Functions

bool finishTrack ()
 finish a track
 
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

StoreArray< MCParticlem_mcParticles
 store array of the MCParticles
 
StoreArray< FANGSSimHitm_simHits
 store array of the SimHits
 
RelationArray m_relMCSimHit {m_mcParticles, m_simHits}
 relation array of the MCParticle -> SimHit relation
 
std::stack< SensorTraversalm_tracks
 Stack of tracks to keep track of particles.
 
Const::EDetector m_subdetector
 Subdetector the class belongs to.
 

Static Private Attributes

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

Sensitive Detector implementation of the FANGS detector.

Definition at line 28 of file SensitiveDetector.h.

Constructor & Destructor Documentation

◆ SensitiveDetector()

Constructor.

Definition at line 30 of file SensitiveDetector.cc.

30 :
31 Simulation::SensitiveDetectorBase("FANGSSensitiveDetector", Const::invalidDetector)
32 {
33 //Register all collections we want to modify and require those we want to use
35 m_simHits.registerInDataStore();
37
38 //Register the Relation so that the TrackIDs get replaced by the actual
39 //MCParticle indices after simulating the events. This is needed as
40 //secondary particles might not be stored so everything relating to those
41 //particles will be attributed to the last saved mother particle
43 }
@ c_negativeWeight
Flip the sign of the weight to become negative if the original element got re-attributed.
Definition: RelationArray.h:79
static void registerMCParticleRelation(const std::string &name, RelationArray::EConsolidationAction ignoreAction=RelationArray::c_negativeWeight)
Register an relation involving MCParticles.
bool registerInDataStore(DataStore::EStoreFlags storeFlags=DataStore::c_WriteOut)
Register the object/array in the DataStore.
RelationArray m_relMCSimHit
relation array of the MCParticle -> SimHit relation
StoreArray< MCParticle > m_mcParticles
store array of the MCParticles
StoreArray< FANGSSimHit > m_simHits
store array of the SimHits

Member Function Documentation

◆ finishTrack()

bool finishTrack ( )
private

finish a track

Definition at line 109 of file SensitiveDetector.cc.

110 {
111 SensorTraversal& traversal = m_tracks.top();
112 // ignore everything below 1eV
113 bool save = traversal.getDepEnergy() > Unit::eV;
114 if (save) {
115 auto momEntry = vecToFloat(traversal.getEntryMomentum());
116 auto posEntry = vecToFloat(traversal.getEntryPosition());
117 auto localposEntry = vecToFloat(traversal.getLocalEntryPosition());
118 auto posExit = vecToFloat(traversal.getExitPosition());
119 int hitIndex = m_simHits.getEntries();
120 m_simHits.appendNew(
121 traversal.getTrackID(),
122 traversal.getLadderID(), traversal.getSensorID(),
123 traversal.getPDGCode(), traversal.getEntryTime(),
124 traversal.getDepEnergy(),
125 traversal.getLength(), posEntry.data(),
126 localposEntry.data(),
127 posExit.data(), momEntry.data()
128 );
129 m_relMCSimHit.add(traversal.getTrackID(), hitIndex, traversal.getDepEnergy());
130 }
131 //No we just need to take care of the stack of traversals
132 if (m_tracks.size() == 1) {
133 //reuse traversal to keep memory if this is the first one
134 traversal.reset();
135 } else {
136 //this should only happen if the parent track got suspended. As this
137 //rarely happens in PXD we do not care for re-usability here
138 m_tracks.pop();
139 }
140 return save;
141 }
void add(index_type from, index_type to, weight_type weight=1.0)
Add a new element to the relation.
static const double eV
[electronvolt]
Definition: Unit.h:112
std::stack< SensorTraversal > m_tracks
Stack of tracks to keep track of particles.

◆ 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 *   
)
overrideprotectedvirtual

Step processing method.

Parameters
stepthe G4Step with the current step information
Returns
true if a Hit has been created, false if the hit was ignored

Implements SensitiveDetectorBase.

Definition at line 45 of file SensitiveDetector.cc.

46 {
47 //Get scintilator ladder and sensor number
48 const G4TouchableHistory* touchable = dynamic_cast<const G4TouchableHistory*>(step->GetPreStepPoint()->GetTouchable());
49 int ladderID = touchable->GetVolume(2)->GetCopyNo();
50 int sensorID = touchable->GetVolume(1)->GetCopyNo();
51
52 //Get Track information
53 const G4Track& track = *step->GetTrack();
54 const int trackID = track.GetTrackID();
55 const int pdgCode = step->GetTrack()->GetDefinition()->GetPDGEncoding();
56 // deposited Energy in Geant4 Units
57 double depEnergy = step->GetTotalEnergyDeposit();
58 // convert into Belle2 Units
59 depEnergy *= Unit::GeV / CLHEP::GeV;
60
61 // get pre and post step point
62 const G4StepPoint& preStep = *step->GetPreStepPoint();
63 const G4StepPoint& postStep = *step->GetPostStepPoint();
64
65 // check if this is the same sensor traversal, otherwise add one to the stack
66 if (m_tracks.empty() || (!m_tracks.top().check(trackID, ladderID, sensorID))) {
67 m_tracks.push(SensorTraversal());
68 }
69 // get the top of the stack
70 SensorTraversal& traversal = m_tracks.top();
71
72 //If new track, remember values
73 if (traversal.getTrackID() == 0) {
74 bool isPrimary = Simulation::TrackInfo::getInfo(track).hasStatus(MCParticle::c_PrimaryParticle);
75 //Add start position
76 const G4ThreeVector preStepPos = preStep.GetPosition() / CLHEP::mm * Unit::mm;
77 const G4ThreeVector preStepMom = preStep.GetMomentum() / CLHEP::MeV * Unit::MeV;
78 //const G4AffineTransform& localToGlobalTransform = preStep.GetTouchableHandle()->GetHistory()->GetTopTransform().Inverse();
79 const G4AffineTransform& localToGlobalTransform = preStep.GetTouchableHandle()->GetHistory()->GetTopTransform();
80 const G4ThreeVector localpreStepPos = localToGlobalTransform.TransformPoint(preStep.GetPosition()) / CLHEP::mm * Unit::mm;
81 const double time = preStep.GetGlobalTime() / CLHEP::ns * Unit::ns;
82 traversal.setInitial(trackID, ladderID, sensorID, pdgCode, isPrimary, preStepPos, localpreStepPos, preStepMom, time);
83 //Remember if the track came from the outside
84 if (preStep.GetStepStatus() == fGeomBoundary) traversal.hasEntered();
85 }
86 //Add current step
87 const G4ThreeVector postStepPos = postStep.GetPosition() / CLHEP::mm * Unit::mm;
88 const double length = step->GetStepLength() / CLHEP::cm * Unit::cm;
89 traversal.add(postStepPos, depEnergy, length);
90
91 //check if we are leaving the volume
92 bool isLeaving = (postStep.GetStepStatus() == fGeomBoundary);
93 if (isLeaving) traversal.hasLeft();
94
95 //if we are leaving or the track is stopped, finish it
96 if (isLeaving || track.GetTrackStatus() >= fStopAndKill) {
97 bool contained = traversal.isContained();
98 bool saved = finishTrack();
99 //we mark all particles as important which created a hit and entered or left the volume
100 if (saved && !contained) {
101 Simulation::TrackInfo::getInfo(track).setIgnore(false);
102 }
103 return saved;
104 }
105 // Track not finished, return false right now
106 return false;
107 }
@ c_PrimaryParticle
bit 0: Particle is primary particle.
Definition: MCParticle.h:47
static const double mm
[millimeters]
Definition: Unit.h:70
static const double MeV
[megaelectronvolt]
Definition: Unit.h:114
static const double cm
Standard units with the value = 1.
Definition: Unit.h:47
static const double ns
Standard of [time].
Definition: Unit.h:48
static const double GeV
Standard of [energy, momentum, mass].
Definition: Unit.h:51
bool step(G4Step *step, G4TouchableHistory *) override
Step processing method.

Member Data Documentation

◆ m_mcParticles

StoreArray<MCParticle> m_mcParticles
private

store array of the MCParticles

Definition at line 43 of file SensitiveDetector.h.

◆ m_relMCSimHit

RelationArray m_relMCSimHit {m_mcParticles, m_simHits}
private

relation array of the MCParticle -> SimHit relation

Definition at line 47 of file SensitiveDetector.h.

◆ m_simHits

StoreArray<FANGSSimHit> m_simHits
private

store array of the SimHits

Definition at line 45 of file SensitiveDetector.h.

◆ m_subdetector

Const::EDetector m_subdetector
privateinherited

Subdetector the class belongs to.

Definition at line 91 of file SensitiveDetectorBase.h.

◆ m_tracks

std::stack<SensorTraversal> m_tracks
private

Stack of tracks to keep track of particles.

Definition at line 50 of file SensitiveDetector.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: