ghit
Loading...
Searching...
No Matches
addHitInfos.cc
Go to the documentation of this file.
1// ghit
2#include "ghit.h"
3
4// glibrary
5#include "gutsConventions.h"
6
7// geant4
8#include "G4VProcess.hh"
9
10// See header for API docs.
11
12void GHit::addHitInfos(const G4Step* step) {
13 auto preStepPoint = step->GetPreStepPoint();
14
15 auto touchable = preStepPoint->GetTouchable();
16
17 // Global position and its local-coordinate transform.
18 G4ThreeVector xyz = preStepPoint->GetPosition();
19 G4ThreeVector xyzL = touchable->GetHistory()->GetTopTransform().TransformPoint(xyz);
20
21 globalPositions.push_back(xyz);
22 localPositions.push_back(xyzL);
23
24 // Energy deposition (scaled by detector multiplier) and global time.
25 double edep = (step->GetTotalEnergyDeposit()) * (gtouchable->getEnergyMultiplier());
26 double time = preStepPoint->GetGlobalTime();
27
28 edeps.push_back(edep);
29 times.push_back(time);
30
31 auto track = step->GetTrack();
32 auto trackVertex = track->GetVertexPosition();
33 int trackId = track->GetTrackID();
34 int motherTrackId = track->GetParentID();
35 int currentPdg = track->GetDefinition()->GetPDGEncoding();
36
37 trackVertexById.emplace(trackId, trackVertex);
38 pdgById.emplace(trackId, currentPdg);
39
40 G4ThreeVector motherTrackVertex(UNINITIALIZEDNUMBERQUANTITY, UNINITIALIZEDNUMBERQUANTITY,
42 int motherPdg = UNINITIALIZEDNUMBERQUANTITY;
43 if (motherTrackId > 0) {
44 auto motherVertex = trackVertexById.find(motherTrackId);
45 if (motherVertex != trackVertexById.end()) {
46 motherTrackVertex = motherVertex->second;
47 }
48 auto motherPdgIt = pdgById.find(motherTrackId);
49 if (motherPdgIt != pdgById.end()) {
50 motherPdg = motherPdgIt->second;
51 }
52 }
53
54 trackVertexPositions.push_back(trackVertex);
55 motherTrackVertexPositions.push_back(motherTrackVertex);
56 pids.push_back(currentPdg);
57 tids.push_back(trackId);
58 motherTids.push_back(motherTrackId);
59 momenta.push_back(preStepPoint->GetMomentum());
60 trackEs.push_back(preStepPoint->GetTotalEnergy());
61 motherPids.push_back(motherPdg);
62
63 if (track->GetCreatorProcess()) {
64 processNames.push_back(track->GetCreatorProcess()->GetProcessName());
65 }
66}
void addHitInfos(const G4Step *thisStep)
Append per-step information from a G4Step.
#define UNINITIALIZEDNUMBERQUANTITY