Belle II Software development
SensitiveDetector.cc
1/**************************************************************************
2 * basf2 (Belle II Analysis Software Framework) *
3 * Author: The Belle II Collaboration *
4 * *
5 * See git log for contributors and copyright holders. *
6 * This file is licensed under LGPL-3.0, see LICENSE.md. *
7 **************************************************************************/
8
9#include <arich/simulation/SensitiveDetector.h>
10#include <arich/dataobjects/ARICHSimHit.h>
11
12// geant4
13#include <G4OpBoundaryProcess.hh>
14#include <G4ProcessManager.hh>
15
16
17#include <simulation/kernel/UserInfo.h>
18#include <G4Track.hh>
19#include <G4Step.hh>
20
21#include <framework/datastore/StoreArray.h>
22#include <framework/datastore/RelationArray.h>
23#include <Math/Vector2D.h>
24#include <TRandom3.h>
25
26using namespace std;
27
28namespace Belle2 {
33 namespace arich {
34
36 Simulation::SensitiveDetectorBase("ARICH", Const::ARICH)
37 {
38
39 StoreArray<MCParticle> particles;
41 RelationArray relation(particles, hits);
43 hits.registerInDataStore();
44 particles.registerRelationTo(hits);
45
46 }
47
48
49 G4bool SensitiveDetector::step(G4Step* aStep, G4TouchableHistory*)
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;
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 }
152
153
154 } // end of namespace arich
156} // end of namespace Belle2
Class ARICHSimHit - Geant4 simulated hit for ARICH.
Definition: ARICHSimHit.h:29
This class provides a set of constants for the framework.
Definition: Const.h:34
Low-level class to create/modify relations between StoreArrays.
Definition: RelationArray.h:62
void add(index_type from, index_type to, weight_type weight=1.0)
Add a new element to the relation.
int getArrayIndex() const
Returns this object's array index (in StoreArray), or -1 if not found.
static void registerMCParticleRelation(const std::string &name, RelationArray::EConsolidationAction ignoreAction=RelationArray::c_negativeWeight)
Register an relation involving MCParticles.
UserInfo class which is used to attach additional information to Geant4 particles and tracks.
Definition: UserInfo.h:36
Accessor to arrays stored in the data store.
Definition: StoreArray.h:113
T * appendNew()
Construct a new T object at the end of the array.
Definition: StoreArray.h:246
DBObjPtr< ARICHSimulationPar > m_simPar
simulation parameters from the DB
bool step(G4Step *aStep, G4TouchableHistory *) override
Process each step and calculate variables defined in PXDSimHit.
Abstract base class for different kinds of events.
STL namespace.