Belle II Software development
SensitiveBar Class Reference

Class providing information on MCParticles hitting the bars. More...

#include <SensitiveBar.h>

Inheritance diagram for SensitiveBar:
SensitiveDetectorBase

Public Member Functions

 SensitiveBar ()
 Constructor.
 
bool step (G4Step *aStep, G4TouchableHistory *) override
 Process each step and fill TOPBarHits.
 
void setReplicaDepth (int depth)
 Sets replica depth of module volume.
 

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.
 

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

int m_replicaDepth = 2
 replica depth of module volume
 
TOPGeometryParm_topgp = TOPGeometryPar::Instance()
 geometry parameters
 
std::vector< int > m_trackIDs
 track ID's
 
StoreArray< MCParticlem_mcParticles
 collection of MC particles
 
StoreArray< TOPBarHitm_barHits
 collection of entrance-to-bar hits
 
RelationArray m_relParticleHit {m_mcParticles, m_barHits}
 relations
 
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

Class providing information on MCParticles hitting the bars.

Applies also quantum efficiency to reduce the number of propagating photons.

Definition at line 31 of file SensitiveBar.h.

Constructor & Destructor Documentation

◆ SensitiveBar()

Constructor.

Definition at line 35 of file SensitiveBar.cc.

35 :
36 Simulation::SensitiveDetectorBase("TOP", Const::TOP)
37 {
38
39 // registration of store arrays and relations
40
41 m_barHits.registerInDataStore();
43
44 // additional registration of MCParticle relation (required for correct relations)
45
47
48 const auto* geo = m_topgp->getGeometry();
49 m_trackIDs.resize(geo->getNumModules(), 0);
50
51 }
@ c_deleteElement
Delete the whole relation element if the original element got re-attributed.
Definition: RelationArray.h:81
static void registerMCParticleRelation(const std::string &name, RelationArray::EConsolidationAction ignoreAction=RelationArray::c_negativeWeight)
Register an relation involving MCParticles.
bool registerRelationTo(const StoreArray< TO > &toArray, DataStore::EDurability durability=DataStore::c_Event, DataStore::EStoreFlags storeFlags=DataStore::c_WriteOut, const std::string &namedRelation="") const
Register a relation to the given StoreArray.
Definition: StoreArray.h:140
std::vector< int > m_trackIDs
track ID's
Definition: SensitiveBar.h:57
RelationArray m_relParticleHit
relations
Definition: SensitiveBar.h:61
TOPGeometryPar * m_topgp
geometry parameters
Definition: SensitiveBar.h:56
StoreArray< MCParticle > m_mcParticles
collection of MC particles
Definition: SensitiveBar.h:59
StoreArray< TOPBarHit > m_barHits
collection of entrance-to-bar hits
Definition: SensitiveBar.h:60
const TOPGeometry * getGeometry() const
Returns pointer to geometry object using basf2 units.

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 }

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

◆ setReplicaDepth()

void setReplicaDepth ( int  depth)
inline

Sets replica depth of module volume.

Parameters
depthreplica depth

Definition at line 51 of file SensitiveBar.h.

51{m_replicaDepth = depth;}
int m_replicaDepth
replica depth of module volume
Definition: SensitiveBar.h:55

◆ step()

bool step ( G4Step *  aStep,
G4TouchableHistory *   
)
overridevirtual

Process each step and fill TOPBarHits.

Parameters
aStepCurrent Geant4 step in the sensitive medium.
Returns
true when particle that is not an optical photon enters the bar

Implements SensitiveDetectorBase.

Definition at line 54 of file SensitiveBar.cc.

55 {
56
57 // get track and particle definition
58
59 G4Track* aTrack = aStep->GetTrack();
60 G4ParticleDefinition* particle = aTrack->GetDefinition();
61
62 // if optical photon, apply QE and return false
63
64 if (particle == G4OpticalPhoton::OpticalPhotonDefinition()) {
65 auto* info = dynamic_cast<Simulation::TrackInfo*>(aTrack->GetUserInformation());
66 if (!info) return false;
67 if (info->getStatus() < 2) {
68 double energy = aTrack->GetKineticEnergy() * Unit::MeV / Unit::eV;
69 double qeffi = m_topgp->getPMTEfficiencyEnvelope(energy);
70 double fraction = info->getFraction();
71 if (qeffi == 0 or gRandom->Uniform() * fraction > qeffi) {
72 aTrack->SetTrackStatus(fStopAndKill);
73 return false;
74 }
75 info->setStatus(2);
76 info->setFraction(qeffi);
77 }
78 return false;
79 }
80
81 // continue for other particles, if they enter the volume for the first time
82
83 G4StepPoint* PrePosition = aStep->GetPreStepPoint();
84 if (PrePosition->GetStepStatus() != fGeomBoundary) return false; // not on boundary
85
86 if (m_barHits.getEntries() == 0) {
87 for (auto& trackID : m_trackIDs) trackID = -1; // reset on new event
88 }
89
90 int moduleID = PrePosition->GetTouchableHandle()->GetReplicaNumber(m_replicaDepth);
91 const auto* geo = m_topgp->getGeometry();
92 if (!geo->isModuleIDValid(moduleID)) {
93 B2ERROR("TOP::SensitiveBar: undefined module ID."
94 << LogVar("moduleID", moduleID));
95 return false;
96 }
97
98 int trackID = aTrack->GetTrackID();
99 if (trackID == m_trackIDs[moduleID - 1]) return false; // not first time on boundary
100 m_trackIDs[moduleID - 1] = trackID;
101
102 // particle other than optical photon entered the volume for the first time
103 // convert to basf2 units and write-out the hit
104
105 G4ThreeVector worldPosition = PrePosition->GetPosition();
106 double tracklength = aTrack->GetTrackLength() - aStep->GetStepLength();
107 double globalTime = PrePosition->GetGlobalTime();
108 G4ThreeVector momentum = PrePosition->GetMomentum();
109
110 ROOT::Math::XYZPoint TPosition(worldPosition.x(), worldPosition.y(), worldPosition.z());
111 ROOT::Math::XYZVector TMomentum(momentum.x(), momentum.y(), momentum.z());
112 ROOT::Math::XYZPoint TOrigin(aTrack->GetVertexPosition().x(),
113 aTrack->GetVertexPosition().y(),
114 aTrack->GetVertexPosition().z());
115
116 // convert to basf2 units
117 TPosition = TPosition * Unit::mm;
118 TMomentum = TMomentum * Unit::MeV;
119 TOrigin = TOrigin * Unit::mm;
120 tracklength = tracklength * Unit::mm;
121
122 const auto& module = geo->getModule(moduleID);
123 auto locPosition = module.pointToLocal(TPosition);
124 auto locMomentum = module.momentumToLocal(TMomentum);
125 double theta = locMomentum.Theta();
126 double phi = locMomentum.Phi();
127
128 int PDG = (int)(particle->GetPDGEncoding());
129
130 // write hit to datastore
131 auto* hit = m_barHits.appendNew(moduleID, PDG, TOrigin, TPosition, TMomentum,
132 globalTime, tracklength, locPosition,
133 theta, phi);
134
135 // set the relation
136 m_relParticleHit.add(trackID, hit->getArrayIndex());
137
138 return true;
139 }
void add(index_type from, index_type to, weight_type weight=1.0)
Add a new element to the relation.
double getPMTEfficiencyEnvelope(double energy) const
Returns PMT efficiency envelope, e.g.
static const double mm
[millimeters]
Definition: Unit.h:70
static const double eV
[electronvolt]
Definition: Unit.h:112
static const double MeV
[megaelectronvolt]
Definition: Unit.h:114
Class to store variables with their name which were sent to the logging service.

Member Data Documentation

◆ m_barHits

StoreArray<TOPBarHit> m_barHits
private

collection of entrance-to-bar hits

Definition at line 60 of file SensitiveBar.h.

◆ m_mcParticles

StoreArray<MCParticle> m_mcParticles
private

collection of MC particles

Definition at line 59 of file SensitiveBar.h.

◆ m_relParticleHit

RelationArray m_relParticleHit {m_mcParticles, m_barHits}
private

relations

Definition at line 61 of file SensitiveBar.h.

◆ m_replicaDepth

int m_replicaDepth = 2
private

replica depth of module volume

Definition at line 55 of file SensitiveBar.h.

◆ m_subdetector

Const::EDetector m_subdetector
privateinherited

Subdetector the class belongs to.

Definition at line 91 of file SensitiveDetectorBase.h.

◆ m_topgp

geometry parameters

Definition at line 56 of file SensitiveBar.h.

◆ m_trackIDs

std::vector<int> m_trackIDs
private

track ID's

Definition at line 57 of file SensitiveBar.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: