Belle II Software development
SensitiveDetector Class Reference

The Class for ARICH Sensitive Detector. More...

#include <SensitiveDetector.h>

Inheritance diagram for SensitiveDetector:
SensitiveDetectorBase

Public Member Functions

 SensitiveDetector ()
 Constructor.
 
bool step (G4Step *aStep, G4TouchableHistory *) override
 Process each step and calculate variables defined in PXDSimHit.
 

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

DBObjPtr< ARICHSimulationParm_simPar
 simulation parameters from the DB
 
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

The Class for ARICH Sensitive Detector.

In this class, every variable defined in ARICHSimHit will be calculated, and stored in datastore.

Definition at line 26 of file SensitiveDetector.h.

Constructor & Destructor Documentation

◆ SensitiveDetector()

Constructor.

Definition at line 35 of file SensitiveDetector.cc.

35 :
36 Simulation::SensitiveDetectorBase("ARICH", Const::ARICH)
37 {
38
39 StoreArray<MCParticle> particles;
40 StoreArray<ARICHSimHit> hits;
41 RelationArray relation(particles, hits);
43 hits.registerInDataStore();
44 particles.registerRelationTo(hits);
45
46 }
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 }

◆ 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()

G4bool step ( G4Step *  aStep,
G4TouchableHistory *   
)
overridevirtual

Process each step and calculate variables defined in PXDSimHit.

Parameters
aStepCurrent Geant4 step in the sensitive medium.
Returns
true if a hit was stored, o.w. false.

Implements SensitiveDetectorBase.

Definition at line 49 of file SensitiveDetector.cc.

50 {
51 //Get particle ID
52 G4Track& track = *aStep->GetTrack();
53 if (track.GetDefinition()->GetParticleName() != "opticalphoton") return false;
54
55 //Get time (check for proper global time)
56 const G4double globalTime = track.GetGlobalTime();
57
58 const G4StepPoint& postStep = *aStep->GetPostStepPoint();
59 const G4ThreeVector& postworldPosition = postStep.GetPosition();
60
61 const G4StepPoint& preStep = *aStep->GetPreStepPoint();
62 const G4ThreeVector& preworldPosition = preStep.GetPosition();
63
64 // direction of photon (+1 is into hapd)
65 int dir = -1;
66 if (postworldPosition.z() - preworldPosition.z() > 0) dir = +1;
67
68 // trigger only on window geometrical boundary
69 if (dir == 1 && postStep.GetStepStatus() != fGeomBoundary) { return false;}
70 if (dir == -1 && preStep.GetStepStatus() != fGeomBoundary) { return false;}
71
72 if ((track.GetNextVolume())->GetName() == "ARICH.HAPDWall") {
73 track.SetTrackStatus(fStopAndKill);
74 return false;
75 }
76
77 // Check if photon is internally reflected in HAPD window
78 G4OpBoundaryProcessStatus theStatus = Undefined;
79
80 G4ProcessManager* OpManager = G4OpticalPhoton::OpticalPhoton()->GetProcessManager();
81
82 if (OpManager) {
83 G4int MAXofPostStepLoops =
84 OpManager->GetPostStepProcessVector()->entries();
85 G4ProcessVector* fPostStepDoItVector =
86 OpManager->GetPostStepProcessVector(typeDoIt);
87 for (G4int i = 0; i < MAXofPostStepLoops; i++) {
88 G4VProcess* fCurrentProcess = (*fPostStepDoItVector)[i];
89 G4OpBoundaryProcess* opProcess = dynamic_cast<G4OpBoundaryProcess*>(fCurrentProcess);
90 if (opProcess) { theStatus = opProcess->GetStatus(); break;}
91 }
92 }
93
94 // if photon is internally reflected and going backward, do nothing
95 if (theStatus == 3 && dir < 0) return 0;
96 if (theStatus == 4 && dir < 0) { if (gRandom->Uniform() < 0.5) track.SetTrackStatus(fStopAndKill); return 0; }
97
98 // apply quantum efficiency if not yet done
99 bool applyQE = true;
100 Simulation::TrackInfo* info =
101 dynamic_cast<Simulation::TrackInfo*>(track.GetUserInformation());
102 if (info) applyQE = info->getStatus() < 2;
103
104 if (applyQE) {
105
106 double energy = track.GetKineticEnergy() / CLHEP::eV;
107 double qeffi = m_simPar->getQE(energy) * m_simPar->getColEff();
108
109 double fraction = 1.;
110 if (info) fraction = info->getFraction();
111 //correct Q.E. for internally reflected photons
112 if (theStatus == 3) qeffi *= m_simPar->getQEScaling();
113 if (gRandom->Uniform() * fraction > qeffi) {
114 // apply possible absorbtion in HAPD window (for internally reflected photons only)
115 if (theStatus == 3 && gRandom->Uniform() < m_simPar->getWindowAbsorbtion()) track.SetTrackStatus(fStopAndKill);
116 return false;
117 }
118 }
119
120 //Get module ID number
121 const G4int moduleID = dir > 0 ? postStep.GetTouchableHandle()->GetReplicaNumber(1) :
122 preStep.GetTouchableHandle()->GetReplicaNumber(2);
123
124 if (moduleID <= 0) return false;
125 //std::cout << "moduleID " << moduleID << std::endl;
126 // std::cout << m_arichgp->getColEff() << " " << info->getFraction() << " " << m_arichgp->getQEScaling() <<" " << m_arichgp->QE(2.0) << std::endl;
127 //Transform to local position
128 G4ThreeVector localPosition = dir > 0 ? postStep.GetTouchableHandle()->GetHistory()->GetTopTransform().TransformPoint(
129 postworldPosition) :
130 preStep.GetTouchableHandle()->GetHistory()->GetTopTransform().TransformPoint(preworldPosition);
131
132 //Get photon energy
133 const G4double energy = track.GetKineticEnergy() / CLHEP::eV;
134
135 //------------------------------------------------------------
136 // Create ARICHSimHit and save it to datastore
137 //------------------------------------------------------------
138
139 ROOT::Math::XYVector locpos(localPosition.x() / CLHEP::cm, localPosition.y() / CLHEP::cm);
140 StoreArray<ARICHSimHit> arichSimHits;
141 ARICHSimHit* simHit = arichSimHits.appendNew(moduleID, locpos, globalTime, energy);
142
143 // add relation to MCParticle
144 StoreArray<MCParticle> mcParticles;
145 RelationArray arichSimHitRel(mcParticles, arichSimHits);
146 arichSimHitRel.add(track.GetParentID(), simHit->getArrayIndex());
147
148 // after detection photon track is killed
149 track.SetTrackStatus(fStopAndKill);
150 return true;
151 }
DBObjPtr< ARICHSimulationPar > m_simPar
simulation parameters from the DB

Member Data Documentation

◆ m_simPar

DBObjPtr<ARICHSimulationPar> m_simPar
private

simulation parameters from the DB

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.

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