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 <beast/srsensor/simulation/SensitiveDetector.h>
10#include <beast/srsensor/dataobjects/SrsensorSimHit.h>
11
12#include <framework/datastore/StoreArray.h>
13#include <framework/datastore/RelationArray.h>
14
15#include <G4Track.hh>
16#include <G4Step.hh>
17
18namespace Belle2 {
24 namespace srsensor {
25
27 Simulation::SensitiveDetectorBase("SrsensorSensitiveDetector", Const::invalidDetector)
28 {
29 //Make sure all collections are registered
30 StoreArray<MCParticle> mcParticles;
32 RelationArray relMCSimHit(mcParticles, simHits);
33
34 //Register all collections we want to modify and require those we want to use
35 mcParticles.registerInDataStore();
36 simHits.registerInDataStore();
37 relMCSimHit.registerInDataStore();
38
39 //Register the Relation so that the TrackIDs get replaced by the actual
40 //MCParticle indices after simulating the events. This is needed as
41 //secondary particles might not be stored so everything relating to those
42 //particles will be attributed to the last saved mother particle
43 registerMCParticleRelation(relMCSimHit);
44 }
45
46 bool SensitiveDetector::step(G4Step* step, G4TouchableHistory*)
47 {
48 //Get Track information
49 const G4Track& track = *step->GetTrack();
50 const int trackID = track.GetTrackID();
51 const double depEnergy = step->GetTotalEnergyDeposit() * CLHEP::MeV;
52 const double nielEnergy = step->GetNonIonizingEnergyDeposit() * CLHEP::MeV;
53 const G4ThreeVector G4tkPos = step->GetTrack()->GetPosition();
54 float tkPos[3];
55 tkPos[0] = G4tkPos.x() * CLHEP::cm;
56 tkPos[1] = G4tkPos.y() * CLHEP::cm;
57 tkPos[2] = G4tkPos.z() * CLHEP::cm;
58 const G4ThreeVector G4tkMom = step->GetTrack()->GetMomentum();
59 float tkMom[3];
60 tkMom[0] = G4tkMom.x() * CLHEP::MeV;
61 tkMom[1] = G4tkMom.y() * CLHEP::MeV;
62 tkMom[2] = G4tkMom.z() * CLHEP::MeV;
63 const G4ThreeVector G4tkMomDir = step->GetTrack()->GetMomentumDirection();
64 float tkMomDir[3];
65 tkMomDir[0] = G4tkMomDir.x() * CLHEP::MeV;
66 tkMomDir[1] = G4tkMomDir.y() * CLHEP::MeV;
67 tkMomDir[2] = G4tkMomDir.z() * CLHEP::MeV;
68 const int tkPDG = step->GetTrack()->GetDefinition()->GetPDGEncoding();
69 const double tkKEnergy = step->GetTrack()->GetKineticEnergy();
70 const int detNb = step->GetTrack()->GetVolume()->GetCopyNo();
71 const double GlTime = step->GetPreStepPoint()->GetGlobalTime();
72 //Ignore everything below 1eV
73 if (depEnergy < CLHEP::eV) return false;
74
75 //Get the datastore arrays
76 StoreArray<MCParticle> mcParticles;
78 RelationArray relMCSimHit(mcParticles, simHits);
79
80 StoreArray<SrsensorSimHit> SrsensorHits;
81 SrsensorSimHit* hit = SrsensorHits.appendNew(
82 depEnergy,
83 nielEnergy,
84 tkPDG,
85 tkKEnergy,
86 detNb,
87 GlTime,
88 tkPos,
89 tkMom,
90 tkMomDir
91 );
92
93 //Add Relation between SimHit and MCParticle with a weight of 1. Since
94 //the MCParticle index is not yet defined we use the trackID from Geant4
95 relMCSimHit.add(trackID, hit->getArrayIndex(), 1.0);
96
97 return true;
98 }
99
100 } //srsensor namespace
102} //Belle2 namespace
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.
static void registerMCParticleRelation(const std::string &name, RelationArray::EConsolidationAction ignoreAction=RelationArray::c_negativeWeight)
Register an relation involving MCParticles.
ClassSrsensorSimHit - Geant4 simulated hit for the Srsensor detector.
bool registerInDataStore(DataStore::EStoreFlags storeFlags=DataStore::c_WriteOut)
Register the object/array in the DataStore.
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
bool step(G4Step *step, G4TouchableHistory *) override
Step processing method.
Abstract base class for different kinds of events.