Belle II Software development
SensitiveDetector Class Reference

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

#include <SensitiveDetector.h>

Inheritance diagram for SensitiveDetector:
SensitiveDetectorBase

Public Member Functions

 SensitiveDetector ()
 Constructor.
 
int saveG4TrackInfo (int trackID, int PDG, float Mass, float Energy, float vtx[3], float mom[3])
 Save saveG4TrackInfo into datastore.
 

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.
 

Protected Attributes

int m_trackID
 track id
 

Private Member Functions

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

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 MICROTPC detector.

Definition at line 22 of file SensitiveDetector.h.

Constructor & Destructor Documentation

◆ SensitiveDetector()

Constructor.

Definition at line 27 of file SensitiveDetector.cc.

27 :
28 Simulation::SensitiveDetectorBase("MicrotpcSensitiveDetector", Const::invalidDetector)
29 {
30 m_trackID = 0;
31 //Make sure all collections are registered
32 StoreArray<MCParticle> mcParticles;
33 StoreArray<MicrotpcSimHit> simHits;
34 StoreArray<TPCG4TrackInfo> TPCG4TrackInfos;
35 RelationArray relMCSimHit(mcParticles, simHits);
36
37 //Register all collections we want to modify and require those we want to use
38 mcParticles.registerInDataStore();
39 simHits.registerInDataStore();
40 relMCSimHit.registerInDataStore();
41 TPCG4TrackInfos.registerInDataStore();
42
43 //Register the Relation so that the TrackIDs get replaced by the actual
44 //MCParticle indices after simulating the events. This is needed as
45 //secondary particles might not be stored so everything relating to those
46 //particles will be attributed to the last saved mother particle
47 registerMCParticleRelation(relMCSimHit);
48 }
static void registerMCParticleRelation(const std::string &name, RelationArray::EConsolidationAction ignoreAction=RelationArray::c_negativeWeight)
Register an relation involving MCParticles.

Member Function Documentation

◆ 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 }

◆ saveG4TrackInfo()

int saveG4TrackInfo ( int  trackID,
int  PDG,
float  Mass,
float  Energy,
float  vtx[3],
float  mom[3] 
)

Save saveG4TrackInfo into datastore.

Definition at line 136 of file SensitiveDetector.cc.

145 {
146 //Get the datastore arrays
147 StoreArray<TPCG4TrackInfo> TPCG4TrackInfos;
148
149 TPCG4TrackInfos.appendNew(TPCG4TrackInfo(trackID, PDG, Mass, Energy, vtx, mom));
150
151 int simhitNumber = TPCG4TrackInfos.getEntries() - 1;
152
153 return (simhitNumber);
154 }//saveG4TrackInfo

◆ 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 50 of file SensitiveDetector.cc.

51 {
52 //Get Track information
53 const G4Track& track = *step->GetTrack();
54 const int trackID = track.GetTrackID();
55 const double depEnergy = step->GetTotalEnergyDeposit() * CLHEP::MeV;
56 const double nielEnergy = step->GetNonIonizingEnergyDeposit() * CLHEP::MeV;
57 const G4ThreeVector G4tkPos = step->GetTrack()->GetPosition();
58 float tkPos[3];
59 tkPos[0] = G4tkPos.x() * CLHEP::cm;
60 tkPos[1] = G4tkPos.y() * CLHEP::cm;
61 tkPos[2] = G4tkPos.z() * CLHEP::cm;
62 const G4ThreeVector G4tkMom = step->GetTrack()->GetMomentum();
63 float tkMom[3];
64 tkMom[0] = G4tkMom.x() * CLHEP::MeV;
65 tkMom[1] = G4tkMom.y() * CLHEP::MeV;
66 tkMom[2] = G4tkMom.z() * CLHEP::MeV;
67 const G4ThreeVector G4tkMomDir = step->GetTrack()->GetMomentumDirection();
68 float tkMomDir[3];
69 tkMomDir[0] = G4tkMomDir.x() * CLHEP::MeV;
70 tkMomDir[1] = G4tkMomDir.y() * CLHEP::MeV;
71 tkMomDir[2] = G4tkMomDir.z() * CLHEP::MeV;
72 const int tkPDG = step->GetTrack()->GetDefinition()->GetPDGEncoding();
73 const double tkKEnergy = step->GetTrack()->GetKineticEnergy();
74 const int detNb = step->GetTrack()->GetVolume()->GetCopyNo();
75 const double GlTime = step->GetPreStepPoint()->GetGlobalTime();
76 //Ignore everything below 1eV
77 //if (depEnergy < CLHEP::eV) return false;
78
79 //Get the datastore arrays
80 StoreArray<MCParticle> mcParticles;
81 StoreArray<MicrotpcSimHit> simHits;
82 RelationArray relMCSimHit(mcParticles, simHits);
83 /*
84 //find out if the process that created the particle was a neutron process
85 bool neuProc = false;
86 G4String CPName;
87 if (step->GetTrack()->GetCreatorProcess() != 0) {
88 const G4VProcess* creator = step->GetTrack()->GetCreatorProcess();
89 CPName = creator->GetProcessName();
90 if (CPName.contains("Neutron")) neuProc = true;
91 }
92 */
93 //set TrackID
94 m_trackID = track.GetTrackID();
95
96 //Save Hit if track leaves volume or is killed
97 if (track.GetNextVolume() != track.GetVolume() || track.GetTrackStatus() >= fStopAndKill) {
98 auto& mcparticle = Simulation::TrackInfo::getInfo(track);
99 int PDG = mcparticle.getPDG();
100 float Mass = mcparticle.getMass();
101 float Energy = mcparticle.getEnergy();
102 float vtx[3];
103 vtx[0] = mcparticle.getProductionVertex().X();
104 vtx[1] = mcparticle.getProductionVertex().Y();
105 vtx[2] = mcparticle.getProductionVertex().Z();
106 float mom[3];
107 mom[0] = mcparticle.getMomentum().X();
108 mom[1] = mcparticle.getMomentum().Y();
109 mom[2] = mcparticle.getMomentum().Z();
110 saveG4TrackInfo(m_trackID, PDG, Mass, Energy, vtx, mom);
111 //Reset TrackID
112 m_trackID = 0;
113 }
114
115 StoreArray<MicrotpcSimHit> MicrotpcHits;
116 MicrotpcSimHit* hit = MicrotpcHits.appendNew(
117 trackID,
118 depEnergy,
119 nielEnergy,
120 tkPDG,
121 tkKEnergy,
122 detNb,
123 GlTime,
124 tkPos,
125 tkMom,
126 tkMomDir
127 );
128
129 //Add Relation between SimHit and MCParticle with a weight of 1. Since
130 //the MCParticle index is not yet defined we use the trackID from Geant4
131 relMCSimHit.add(trackID, hit->getArrayIndex(), 1.0);
132
133 return true;
134 }
int saveG4TrackInfo(int trackID, int PDG, float Mass, float Energy, float vtx[3], float mom[3])
Save saveG4TrackInfo into datastore.
bool step(G4Step *step, G4TouchableHistory *) override
Step processing method.

Member Data Documentation

◆ m_subdetector

Const::EDetector m_subdetector
privateinherited

Subdetector the class belongs to.

Definition at line 91 of file SensitiveDetectorBase.h.

◆ m_trackID

int m_trackID
protected

track id

Definition at line 45 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: