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/plume/simulation/SensitiveDetector.h>
10#include <beast/plume/dataobjects/PlumeSimHit.h>
11#include <beast/plume/dataobjects/PlumeHit.h>
12
13#include <framework/datastore/StoreArray.h>
14#include <framework/datastore/RelationArray.h>
15
16#include <G4Track.hh>
17#include <G4Step.hh>
18
19#include <iostream>
20#define _USE_MATH_DEFINES
21#include <cmath>
22using namespace std;
23
24
25namespace Belle2 {
31 namespace plume {
32
34 Simulation::SensitiveDetectorBase("PlumeSensitiveDetector", Const::invalidDetector)
35 {
40 current_pdgID = 0;
42 current_posIN_u = 0.;
43 current_posIN_v = 0.;
44 current_posIN_w = 0.;
45 current_posIN_x = 0.;
46 current_posIN_y = 0.;
47 current_posIN_z = 0.;
58 current_nielDep = 0.;
59 current_thetaAngle = 0.; // local sensor frame
60 current_phiAngle = 0.; // local sensor frame
62
63 //Make sure all collections are registered
64 StoreArray<MCParticle> mcParticles;
66 StoreArray<PlumeHit> plumeHits;
67
68 RelationArray relMCSimHit(mcParticles, simHits);
69 RelationArray relMCplumeHit(mcParticles, plumeHits);
70
71 //Register all collections we want to modify and require those we want to use
72 mcParticles.registerInDataStore();
73 simHits.registerInDataStore();
74 plumeHits.registerInDataStore();
75 relMCSimHit.registerInDataStore();
76 relMCplumeHit.registerInDataStore();
77
78 //Register the Relation so that the TrackIDs get replaced by the actual
79 //MCParticle indices after simulating the events. This is needed as
80 //secondary particles might not be stored so everything relating to those
81 //particles will be attributed to the last saved mother particle
82 registerMCParticleRelation(relMCSimHit);
83 registerMCParticleRelation(relMCplumeHit);
84 }
85
87 {
88
89 }
90
91
92 bool SensitiveDetector::step(G4Step* step, G4TouchableHistory*)
93 {
94
95 //Get track information, and pre- and post-step
96 const G4Track& track = *step->GetTrack();
97 const G4StepPoint& preStepPoint = *step->GetPreStepPoint();
98 const G4StepPoint& postStepPoint = *step->GetPostStepPoint();
99
100
101 // If new track, store general information from this first step
102 if (current_trackID != track.GetTrackID()) { // if new track
103 current_trackID = track.GetTrackID();
104 current_trackVertex_x = track.GetVertexPosition().x();
105 current_trackVertex_y = track.GetVertexPosition().y();
106 current_trackVertex_z = track.GetVertexPosition().z();
107 current_pdgID = track.GetDefinition()->GetPDGEncoding();
108 current_sensorID = preStepPoint.GetTouchableHandle()->GetReplicaNumber(1);
109
110 // since this is first step in volume store track incidence and momentum and reset energy loss
111 G4ThreeVector preStepPointPosition = preStepPoint.GetPosition();
112 G4ThreeVector trackMomentum = track.GetMomentum();
114 current_nielDep = 0.;
115
116 current_posIN_x = preStepPointPosition.x();
117 current_posIN_y = preStepPointPosition.y();
118 current_posIN_z = preStepPointPosition.z();
119 current_momentum_x = trackMomentum.x();
120 current_momentum_y = trackMomentum.y();
121 current_momentum_z = trackMomentum.z();
122
123 // We want the track angles of incidence on the local volume and the local position
124 G4TouchableHandle theTouchable = preStepPoint.GetTouchableHandle();
125 G4ThreeVector worldDirection = preStepPoint.GetMomentumDirection();
126 G4ThreeVector localDirection = theTouchable->GetHistory()->GetTopTransform().TransformAxis(worldDirection);
127
128 if (localDirection.z() < 0.) {
129 current_thetaAngle = M_PI - acos(-1. * localDirection.z());
130 } else {
131 current_thetaAngle = acos(localDirection.z());
132 }
133
134
135 if (localDirection.x() < 0. && localDirection.y() > 0.) {
136 current_phiAngle = atan(-1. * localDirection.x() / localDirection.y());
137 } else if (localDirection.x() < 0. && localDirection.y() < 0.) {
138 current_phiAngle = M_PI - atan(localDirection.x() / localDirection.y());
139 } else if (localDirection.x() > 0. && localDirection.y() < 0.) {
140 current_phiAngle = M_PI + atan(-1. * localDirection.x() / localDirection.y());
141 } else if (localDirection.x() > 0. && localDirection.y() > 0.) {
142 current_phiAngle = 2. * M_PI - atan(localDirection.x() / localDirection.y());
143 } else if (localDirection.x() < 0. && localDirection.y() == 0.) {
144 current_phiAngle = M_PI / 2.;
145 } else if (localDirection.x() > 0. && localDirection.y() == 0.) {
146 current_phiAngle = 3.* M_PI / 2.;
147 } else if (localDirection.x() == 0. && localDirection.y() > 0.) {
148 current_phiAngle = 0.;
149 } else if (localDirection.x() == 0. && localDirection.y() < 0.) {
150 current_phiAngle = M_PI;
151 }
152
153 G4ThreeVector localINPosition = theTouchable->GetHistory()->GetTopTransform().TransformPoint(preStepPointPosition);
154 current_posIN_u = localINPosition.y();
155 current_posIN_v = - localINPosition.x();
156 current_posIN_w = localINPosition.z();
157
158 current_globalTime = step->GetPreStepPoint()->GetGlobalTime();
159
160 } // end if new track
161
162 // Update information at each step
163 current_energyDep += step->GetTotalEnergyDeposit();
164 current_nielDep += step->GetNonIonizingEnergyDeposit();
165
166 // If track leaves volume or is killed, store final step info and save simHit
167 if (track.GetNextVolume() != track.GetVolume() || track.GetTrackStatus() >= fStopAndKill) { // if last step
168
169 G4ThreeVector postStepPointPosition = postStepPoint.GetPosition();
170
171 current_posOUT_x = postStepPointPosition.x();
172 current_posOUT_y = postStepPointPosition.y();
173 current_posOUT_z = postStepPointPosition.z();
174 G4ThreeVector localOUTPosition = preStepPoint.GetTouchableHandle()->GetHistory()->GetTopTransform().TransformPoint(
175 postStepPointPosition);
176 current_posOUT_u = localOUTPosition.y();
177 current_posOUT_v = - localOUTPosition.x();
178 current_posOUT_w = localOUTPosition.z();
179
180 //Get the datastore arrays
181 StoreArray<MCParticle> mcParticles;
183
184 if (current_energyDep > CLHEP::eV) {
185
186 RelationArray relMCSimHit(mcParticles, simHits);
187 PlumeSimHit* hit = simHits.appendNew(
191 current_trackVertex_x / CLHEP::mm,
192 current_trackVertex_y / CLHEP::mm,
193 current_trackVertex_z / CLHEP::mm,
194 current_energyDep / CLHEP::MeV,
195 current_nielDep / CLHEP::MeV,
196 current_posIN_x / CLHEP::mm,
197 current_posIN_y / CLHEP::mm,
198 current_posIN_z / CLHEP::mm,
199 current_posIN_u / CLHEP::mm,
200 current_posIN_v / CLHEP::mm,
201 current_posIN_w / CLHEP::mm,
202 current_posOUT_u / CLHEP::mm,
203 current_posOUT_v / CLHEP::mm,
204 current_posOUT_w / CLHEP::mm,
205 current_momentum_x / CLHEP::GeV,
206 current_momentum_y / CLHEP::GeV,
207 current_momentum_z / CLHEP::GeV,
208 current_thetaAngle / CLHEP::degree,
209 current_phiAngle / CLHEP::degree,
210 current_globalTime / CLHEP::nanosecond
211 );
212
213 //Add Relation between SimHit and MCParticle with a weight of 1. Since
214 //the MCParticle index is not yet defined we use the trackID from Geant4
215 relMCSimHit.add(current_trackID, hit->getArrayIndex(), 1.0);
216 }
217
218 //Reset TrackID
219 current_trackID = 0;
220
221 } // end if last
222
223
224 //Ignore everything below 1eV
225 if (current_energyDep < CLHEP::eV) return false;
226
227
228 return true;
229 }
230
231 } //plume namespace
233} //Belle2 namespace
This class provides a set of constants for the framework.
Definition: Const.h:34
Class PlumeSimHit - Geant4 simulated hit for the PLUME detector.
Definition: PlumeSimHit.h:27
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.
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
float current_trackVertex_z
track production vertex z coordinates in G4 ref
float current_posIN_v
incoming track position v in sensor ref, in mm
float current_posIN_u
incoming track position u in sensor ref, in mm
float current_momentum_z
incoming track momentum, z coordinates in G4 ref, in GeV
float current_posIN_y
incoming track position y in G4 ref, in mm
float current_momentum_y
incoming track momentum, y coordinates in G4 ref, in GeV
float current_posOUT_w
outcoming track position w in sensor ref, in mm
float current_energyDep
deposited energy in MeV
float current_trackVertex_x
track production vertex x coordinates in G4 ref
float current_posOUT_z
outcoming track position z in sensor ref, in mm
float current_posOUT_v
outcoming track position v in sensor ref, in mm
float current_posIN_z
incoming track position z in G4 ref, in mm
float current_thetaAngle
local (sensor ref) theta angle, out of sensor plane, in degree
float current_posOUT_u
outcoming track position u in sensor ref, in mm
float current_momentum_x
incoming track momentum, x coordinates in G4 ref, in GeV
float current_trackVertex_y
track production vertex y coordinates in G4 ref
float current_posOUT_y
outcoming track position y in sensor ref, in mm
float current_posIN_x
incoming track position x in G4 ref, in mm
float current_phiAngle
local (sensor ref) phi angle, in sensor plane, in degree
float current_nielDep
non ionizing deposited energy
float current_posIN_w
incoming track position w in sensor ref, in mm
bool step(G4Step *step, G4TouchableHistory *) override
Step processing method.
float current_posOUT_x
outcoming track position x in sensor ref, in mm
double atan(double a)
atan for double
Definition: beamHelpers.h:34
Abstract base class for different kinds of events.
STL namespace.