Belle II Software  release-06-01-15
SensitiveDetector Class Reference

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

#include <SensitiveDetector.h>

Inheritance diagram for SensitiveDetector:
Collaboration diagram for SensitiveDetector:

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. More...
 
static void registerMCParticleRelation (const std::string &name, RelationArray::EConsolidationAction ignoreAction=RelationArray::c_negativeWeight)
 Register an relation involving MCParticles. More...
 
static void registerMCParticleRelation (const RelationArray &relation, RelationArray::EConsolidationAction ignoreAction=RelationArray::c_negativeWeight)
 Overload to make it easer to register MCParticle relations. More...
 

Protected Member Functions

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

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. More...
 

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.

Member Function Documentation

◆ 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 96 of file SensitiveDetectorBase.h.

◆ 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 68 of file SensitiveDetectorBase.h.

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

◆ 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 52 of file SensitiveDetectorBase.h.

◆ 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  if (m_trackID != track.GetTrackID()) {
94  //TrackID changed, store track informations
95  m_trackID = track.GetTrackID();
96  }
97 
98  //Save Hit if track leaves volume or is killed
99  if (track.GetNextVolume() != track.GetVolume() || track.GetTrackStatus() >= fStopAndKill) {
100  auto& mcparticle = Simulation::TrackInfo::getInfo(track);
101  int PDG = mcparticle.getPDG();
102  float Mass = mcparticle.getMass();
103  float Energy = mcparticle.getEnergy();
104  float vtx[3];
105  vtx[0] = mcparticle.getProductionVertex().X();
106  vtx[1] = mcparticle.getProductionVertex().Y();
107  vtx[2] = mcparticle.getProductionVertex().Z();
108  float mom[3];
109  mom[0] = mcparticle.getMomentum().X();
110  mom[1] = mcparticle.getMomentum().Y();
111  mom[2] = mcparticle.getMomentum().Z();
112  saveG4TrackInfo(m_trackID, PDG, Mass, Energy, vtx, mom);
113  //Reset TrackID
114  m_trackID = 0;
115  }
116 
117  StoreArray<MicrotpcSimHit> MicrotpcHits;
118  MicrotpcSimHit* hit = MicrotpcHits.appendNew(
119  trackID,
120  depEnergy,
121  nielEnergy,
122  tkPDG,
123  tkKEnergy,
124  detNb,
125  GlTime,
126  tkPos,
127  tkMom,
128  tkMomDir
129  );
130 
131  //Add Relation between SimHit and MCParticle with a weight of 1. Since
132  //the MCParticle index is not yet defined we use the trackID from Geant4
133  relMCSimHit.add(trackID, hit->getArrayIndex(), 1.0);
134 
135  return true;
136  }
static Payload getInfo(Carrier &obj)
Static function to just return UserInformation attached to the obj of type Carrier.
Definition: UserInfo.h:100
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.

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