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/microtpc/simulation/SensitiveDetector.h>
10#include <beast/microtpc/dataobjects/MicrotpcSimHit.h>
11#include <beast/microtpc/dataobjects/TPCG4TrackInfo.h>
12
13#include <framework/datastore/StoreArray.h>
14#include <framework/datastore/RelationArray.h>
15
16#include <G4Track.hh>
17#include <G4Step.hh>
18
19namespace Belle2 {
25 namespace microtpc {
26
28 Simulation::SensitiveDetectorBase("MicrotpcSensitiveDetector", Const::invalidDetector)
29 {
30 m_trackID = 0;
31 //Make sure all collections are registered
32 StoreArray<MCParticle> mcParticles;
34 StoreArray<TPCG4TrackInfo> TPCG4TrackInfos;
35 RelationArray relMCSimHit(mcParticles, simHits);
36
37 //Register all collections we want to modify and require those we want to use
38 mcParticles.registerInDataStore();
39 simHits.registerInDataStore();
40 relMCSimHit.registerInDataStore();
41 TPCG4TrackInfos.registerInDataStore();
42
43 //Register the Relation so that the TrackIDs get replaced by the actual
44 //MCParticle indices after simulating the events. This is needed as
45 //secondary particles might not be stored so everything relating to those
46 //particles will be attributed to the last saved mother particle
47 registerMCParticleRelation(relMCSimHit);
48 }
49
50 bool SensitiveDetector::step(G4Step* step, G4TouchableHistory*)
51 {
52 //Get Track information
53 const G4Track& track = *step->GetTrack();
54 const int trackID = track.GetTrackID();
55 const double depEnergy = step->GetTotalEnergyDeposit() * CLHEP::MeV;
56 const double nielEnergy = step->GetNonIonizingEnergyDeposit() * CLHEP::MeV;
57 const G4ThreeVector G4tkPos = step->GetTrack()->GetPosition();
58 float tkPos[3];
59 tkPos[0] = G4tkPos.x() * CLHEP::cm;
60 tkPos[1] = G4tkPos.y() * CLHEP::cm;
61 tkPos[2] = G4tkPos.z() * CLHEP::cm;
62 const G4ThreeVector G4tkMom = step->GetTrack()->GetMomentum();
63 float tkMom[3];
64 tkMom[0] = G4tkMom.x() * CLHEP::MeV;
65 tkMom[1] = G4tkMom.y() * CLHEP::MeV;
66 tkMom[2] = G4tkMom.z() * CLHEP::MeV;
67 const G4ThreeVector G4tkMomDir = step->GetTrack()->GetMomentumDirection();
68 float tkMomDir[3];
69 tkMomDir[0] = G4tkMomDir.x() * CLHEP::MeV;
70 tkMomDir[1] = G4tkMomDir.y() * CLHEP::MeV;
71 tkMomDir[2] = G4tkMomDir.z() * CLHEP::MeV;
72 const int tkPDG = step->GetTrack()->GetDefinition()->GetPDGEncoding();
73 const double tkKEnergy = step->GetTrack()->GetKineticEnergy();
74 const int detNb = step->GetTrack()->GetVolume()->GetCopyNo();
75 const double GlTime = step->GetPreStepPoint()->GetGlobalTime();
76 //Ignore everything below 1eV
77 //if (depEnergy < CLHEP::eV) return false;
78
79 //Get the datastore arrays
80 StoreArray<MCParticle> mcParticles;
82 RelationArray relMCSimHit(mcParticles, simHits);
83 /*
84 //find out if the process that created the particle was a neutron process
85 bool neuProc = false;
86 G4String CPName;
87 if (step->GetTrack()->GetCreatorProcess() != 0) {
88 const G4VProcess* creator = step->GetTrack()->GetCreatorProcess();
89 CPName = creator->GetProcessName();
90 if (CPName.contains("Neutron")) neuProc = true;
91 }
92 */
93 //set TrackID
94 m_trackID = track.GetTrackID();
95
96 //Save Hit if track leaves volume or is killed
97 if (track.GetNextVolume() != track.GetVolume() || track.GetTrackStatus() >= fStopAndKill) {
98 auto& mcparticle = Simulation::TrackInfo::getInfo(track);
99 int PDG = mcparticle.getPDG();
100 float Mass = mcparticle.getMass();
101 float Energy = mcparticle.getEnergy();
102 float vtx[3];
103 vtx[0] = mcparticle.getProductionVertex().X();
104 vtx[1] = mcparticle.getProductionVertex().Y();
105 vtx[2] = mcparticle.getProductionVertex().Z();
106 float mom[3];
107 mom[0] = mcparticle.getMomentum().X();
108 mom[1] = mcparticle.getMomentum().Y();
109 mom[2] = mcparticle.getMomentum().Z();
110 saveG4TrackInfo(m_trackID, PDG, Mass, Energy, vtx, mom);
111 //Reset TrackID
112 m_trackID = 0;
113 }
114
115 StoreArray<MicrotpcSimHit> MicrotpcHits;
116 MicrotpcSimHit* hit = MicrotpcHits.appendNew(
117 trackID,
118 depEnergy,
119 nielEnergy,
120 tkPDG,
121 tkKEnergy,
122 detNb,
123 GlTime,
124 tkPos,
125 tkMom,
126 tkMomDir
127 );
128
129 //Add Relation between SimHit and MCParticle with a weight of 1. Since
130 //the MCParticle index is not yet defined we use the trackID from Geant4
131 relMCSimHit.add(trackID, hit->getArrayIndex(), 1.0);
132
133 return true;
134 }
135
137 (
138 int trackID,
139 int PDG,
140 float Mass,
141 float Energy,
142 float vtx[3],
143 float mom[3]
144 )
145 {
146 //Get the datastore arrays
147 StoreArray<TPCG4TrackInfo> TPCG4TrackInfos;
148
149 TPCG4TrackInfos.appendNew(TPCG4TrackInfo(trackID, PDG, Mass, Energy, vtx, mom));
150
151 int simhitNumber = TPCG4TrackInfos.getEntries() - 1;
152
153 return (simhitNumber);
154 }//saveG4TrackInfo
155
156 } //microtpc namespace
158} //Belle2 namespace
This class provides a set of constants for the framework.
Definition: Const.h:34
ClassMicrotpcSimHit - Geant4 simulated hit for the Microtpc detector.
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.
static Payload getInfo(Carrier &obj)
Static function to just return UserInformation attached to the obj of type Carrier.
Definition: UserInfo.h:100
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
int getEntries() const
Get the number of objects in the array.
Definition: StoreArray.h:216
Class TPCG4TrackInfo - Geant4 simulated hit for the Microtpc detector.
int saveG4TrackInfo(int trackID, int PDG, float Mass, float Energy, float vtx[3], float mom[3])
Save saveG4TrackInfo into datastore.
bool step(G4Step *step, G4TouchableHistory *) override
Step processing method.
Abstract base class for different kinds of events.