Belle II Software  release-06-02-00
SensitivePMT.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 <top/simulation/SensitivePMT.h>
10 
11 #include <simulation/kernel/UserInfo.h>
12 #include <G4Step.hh>
13 #include <G4Track.hh>
14 #include <G4UnitsTable.hh>
15 #include <G4ParticleDefinition.hh>
16 #include <G4ParticleTypes.hh>
17 
18 #include <framework/logging/Logger.h>
19 #include <framework/gearbox/Unit.h>
20 
21 #include <TVector3.h>
22 #include <TRandom.h>
23 
24 using namespace std;
25 
26 namespace Belle2 {
31  namespace TOP {
32 
33  SensitivePMT::SensitivePMT():
34  Simulation::SensitiveDetectorBase("TOP", Const::TOP)
35  {
36 
37  m_simHits.registerInDataStore();
39 
40  m_simPhotons.registerInDataStore(DataStore::c_DontWriteOut);
41  m_simHits.registerRelationTo(m_simPhotons, DataStore::c_Event,
43 
45 
46  }
47 
48 
49  G4bool SensitivePMT::step(G4Step* aStep, G4TouchableHistory*)
50  {
51  // photon track
52  G4Track& photon = *aStep->GetTrack();
53 
54  // check if the track is an optical photon
55  if (photon.GetDefinition() != G4OpticalPhoton::OpticalPhotonDefinition()) return false;
56 
57  // photon energy in [eV]
58  double energy = photon.GetKineticEnergy() * Unit::MeV / Unit::eV;
59 
60  // pmt and module ID
61  int pmtID = photon.GetTouchableHandle()->GetReplicaNumber(m_pmtReplicaDepth);
62  int moduleID = photon.GetTouchableHandle()->GetReplicaNumber(m_moduleReplicaDepth);
63 
64  // hit position in local frame, converted to Basf units
65  G4ThreeVector localPosition = photon.GetTouchableHandle()->GetHistory()->GetTopTransform().TransformPoint(photon.GetPosition());
66  double xLocal = localPosition.x() * Unit::mm;
67  double yLocal = localPosition.y() * Unit::mm;
68 
69  // apply quantum efficiency if not yet done
70  bool applyQE = true;
71  double fraction = 1;
72  auto* info = dynamic_cast<Simulation::TrackInfo*>(photon.GetUserInformation());
73  if (info) {
74  applyQE = info->getStatus() < 3;
75  fraction = info->getFraction();
76  }
77  if (applyQE) {
78  double qeffi = m_topgp->getPMTEfficiency(energy, moduleID, pmtID, xLocal, yLocal);
79  if (qeffi == 0 or gRandom->Uniform() * fraction > qeffi) {
80  photon.SetTrackStatus(fStopAndKill);
81  return false;
82  }
83  info->setStatus(3);
84  info->setFraction(qeffi);
85  }
86 
87  // photon at detection
88  const G4ThreeVector& g_detPoint = photon.GetPosition();
89  const G4ThreeVector& g_detMomDir = photon.GetMomentumDirection();
90  TVector3 detPoint(g_detPoint.x(), g_detPoint.y(), g_detPoint.z());
91  TVector3 detMomDir(g_detMomDir.x(), g_detMomDir.y(), g_detMomDir.z());
92  double detTime = photon.GetGlobalTime();
93  double length = photon.GetTrackLength();
94 
95  // photon at emission
96  const G4ThreeVector& g_emiPoint = photon.GetVertexPosition();
97  const G4ThreeVector& g_emiMomDir = photon.GetVertexMomentumDirection();
98  TVector3 emiPoint(g_emiPoint.x(), g_emiPoint.y(), g_emiPoint.z());
99  TVector3 emiMomDir(g_emiMomDir.x(), g_emiMomDir.y(), g_emiMomDir.z());
100  double emiTime = photon.GetGlobalTime() - photon.GetLocalTime();
101 
102  // convert to basf2 units
103  emiPoint = emiPoint * Unit::mm;
104  detPoint = detPoint * Unit::mm;
105  length = length * Unit::mm;
106 
107  // write to store arrays; add relations
108 
109  auto* simHit = m_simHits.appendNew(moduleID, pmtID, xLocal, yLocal,
110  detTime, energy);
111 
112  int parentID = photon.GetParentID();
113  if (parentID == 0) parentID = photon.GetTrackID();
114  m_relParticleHit.add(parentID, simHit->getArrayIndex());
115 
116  const auto* geo = m_topgp->getGeometry();
117  if (geo->isModuleIDValid(moduleID)) {
118  // transform to local frame
119  const auto& module = geo->getModule(moduleID);
120  emiPoint = module.pointToLocal(emiPoint);
121  detPoint = module.pointToLocal(detPoint);
122  emiMomDir = module.momentumToLocal(emiMomDir);
123  detMomDir = module.momentumToLocal(detMomDir);
124  } else {
125  B2ERROR("TOP::SensitivePMT: undefined module ID."
126  << LogVar("moduleID", moduleID));
127  }
128  auto* simPhoton = m_simPhotons.appendNew(moduleID,
129  emiPoint, emiMomDir, emiTime,
130  detPoint, detMomDir, detTime,
131  length, energy);
132 
133  simHit->addRelationTo(simPhoton);
134 
135  // kill photon after detection
136  photon.SetTrackStatus(fStopAndKill);
137 
138  return true;
139  }
140 
141 
142  } // end of namespace top
144 } // end of namespace Belle2
This class provides a set of constants for the framework.
Definition: Const.h:34
@ c_DontWriteOut
Object/array should be NOT saved by output modules.
Definition: DataStore.h:71
@ c_Event
Different object in each event, all objects/arrays are invalidated after event() function has been ca...
Definition: DataStore.h:59
void add(index_type from, index_type to, weight_type weight=1.0)
Add a new element to the relation.
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
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
RelationArray m_relParticleHit
relations
Definition: SensitivePMT.h:67
int m_pmtReplicaDepth
replica depth of PMT volume
Definition: SensitivePMT.h:60
TOPGeometryPar * m_topgp
geometry parameters
Definition: SensitivePMT.h:62
int m_moduleReplicaDepth
replica depth of module volume
Definition: SensitivePMT.h:61
StoreArray< MCParticle > m_mcParticles
collection of MC particles
Definition: SensitivePMT.h:64
bool step(G4Step *aStep, G4TouchableHistory *) override
Process each step, fill TOPSimHits and TOPSimPhotons.
Definition: SensitivePMT.cc:49
StoreArray< TOPSimPhoton > m_simPhotons
collection of simulated photons
Definition: SensitivePMT.h:66
StoreArray< TOPSimHit > m_simHits
collection of simulated hits
Definition: SensitivePMT.h:65
const TOPGeometry * getGeometry() const
Returns pointer to geometry object using basf2 units.
double getPMTEfficiency(double energy, int moduleID, int pmtID, double x, double y) const
Returns PMT pixel efficiency, a product of quantum and collection efficiency.
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.
Abstract base class for different kinds of events.