Belle II Software  release-06-02-00
EKLMSensitiveDetector.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 /* Own header. */
10 #include <klm/eklm/simulation/EKLMSensitiveDetector.h>
11 
12 /* KLM headers. */
13 #include <klm/dbobjects/eklm/EKLMSimulationParameters.h>
14 
15 /* Belle 2 headers. */
16 #include <framework/database/DBObjPtr.h>
17 #include <framework/gearbox/Unit.h>
18 #include <framework/logging/Logger.h>
19 
20 /* Geant4 headers. */
21 #include <G4Step.hh>
22 
23 /* CLHEP headers. */
24 #include <CLHEP/Geometry/Point3D.h>
25 
26 using namespace Belle2;
27 
29  Simulation::SensitiveDetectorBase(name, Const::KLM),
30  m_ElementNumbers(&(EKLMElementNumbers::Instance()))
31 {
33  if (!simPar.isValid())
34  B2FATAL("EKLM simulation parameters are not available.");
35  m_ThresholdHitTime = simPar->getHitTimeThreshold();
37  m_SimHits.registerInDataStore();
40 }
41 
43 {
44 }
45 
46 bool EKLM::EKLMSensitiveDetector::step(G4Step* aStep, G4TouchableHistory*)
47 {
48  const int stripLevel = 1;
49  int section, layer, sector, plane, strip, stripGlobal;
50  HepGeom::Point3D<double> gpos, lpos;
51  G4TouchableHandle hist = aStep->GetPreStepPoint()->GetTouchableHandle();
52  section = hist->GetVolume(stripLevel + 6)->GetCopyNo();
53  layer = hist->GetVolume(stripLevel + 5)->GetCopyNo();
54  sector = hist->GetVolume(stripLevel + 4)->GetCopyNo();
55  plane = hist->GetVolume(stripLevel + 3)->GetCopyNo();
56  strip = hist->GetVolume(stripLevel)->GetCopyNo();
57  stripGlobal = m_ElementNumbers->stripNumber(
58  section, layer, sector, plane, strip);
59  const G4double eDep = aStep->GetTotalEnergyDeposit();
60  /* Do not record hits without deposited energy. */
61  if (eDep <= 0)
62  return false;
63  const G4Track& track = * aStep->GetTrack();
64  const G4double hitTime = track.GetGlobalTime();
65  /* No time cut for background studies. */
66  if (hitTime > m_ThresholdHitTime) {
67  B2INFO("EKLMSensitiveDetector: "
68  " ALL HITS WITH TIME > hitTimeThreshold ARE DROPPED!!");
69  return false;
70  }
71  /* Hit position. */
72  gpos = 0.5 * (aStep->GetPostStepPoint()->GetPosition() +
73  aStep->GetPreStepPoint()->GetPosition());
74  lpos = hist->GetHistory()->GetTopTransform().TransformPoint(gpos);
75  /* Create step hit and store in to DataStore */
76  EKLMSimHit* hit = m_SimHits.appendNew();
77  CLHEP::Hep3Vector trackMomentum = track.GetMomentum();
78  hit->setMomentum(TLorentzVector(trackMomentum.x(), trackMomentum.y(),
79  trackMomentum.z(), track.GetTotalEnergy()));
80  hit->setTrackID(track.GetTrackID());
81  hit->setParentTrackID(track.GetParentID());
82  hit->setLocalPosition(lpos.x() / CLHEP::mm * Unit::mm,
83  lpos.y() / CLHEP::mm * Unit::mm,
84  lpos.z() / CLHEP::mm * Unit::mm);
85  hit->setPosition(gpos.x() / CLHEP::mm * Unit::mm,
86  gpos.y() / CLHEP::mm * Unit::mm,
87  gpos.z() / CLHEP::mm * Unit::mm);
88  hit->setEnergyDeposit(eDep);
89  hit->setPDG(track.GetDefinition()->GetPDGEncoding());
90  hit->setTime(hitTime);
91  hit->setStrip(strip);
92  hit->setPlane(plane);
93  hit->setSector(sector);
94  hit->setLayer(layer);
95  hit->setSection(section);
96  hit->setVolumeID(stripGlobal);
97  m_MCParticlesToSimHits.add(track.GetTrackID(), hit->getArrayIndex());
98  return true;
99 }
100 
This class provides a set of constants for the framework.
Definition: Const.h:34
bool isValid() const
Check whether a valid object was obtained from the database.
Class for accessing objects in the database.
Definition: DBObjPtr.h:21
EKLM element numbers.
Class EKLMSimHit stores information on particular Geant step; using information from TrackID and Pare...
Definition: EKLMSimHit.h:32
bool step(G4Step *aStep, G4TouchableHistory *history) override
Process each step and calculate variables for EKLMSimHit store EKLMSimHit.
EKLMSensitiveDetector(G4String name)
Constructor.
StoreArray< EKLMSimHit > m_SimHits
Simulation hits.
G4double m_ThresholdHitTime
All hits with time large than m_ThresholdHitTime will be dropped.
RelationArray m_MCParticlesToSimHits
Relation array between MCPartices and EKLMSimHits.
StoreArray< MCParticle > m_MCParticles
MC particles.
static void registerMCParticleRelation(const std::string &name, RelationArray::EConsolidationAction ignoreAction=RelationArray::c_negativeWeight)
Register an relation involving MCParticles.
bool isOptional(const std::string &name="")
Tell the DataStore about an optional input.
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
static const double mm
[millimeters]
Definition: Unit.h:70
Abstract base class for different kinds of events.