ghit
Loading...
Searching...
No Matches
ghit.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 "G4VVisManager.hh"
9#include "G4Circle.hh"
10#include "G4VisAttributes.hh"
11#include "Randomize.hh"
12
13using std::string;
14using std::vector;
15
16std::atomic<int> GHit::globalHitCounter{0};
17thread_local std::map<int, G4ThreeVector> GHit::trackVertexById;
18thread_local std::map<int, int> GHit::pdgById;
19
20// MT definitions, as from:
21// https://twiki.cern.ch/twiki/bin/view/Geant4/QuickMigrationGuideForGeant4V10
22G4ThreadLocal G4Allocator<GHit>* GHitAllocator = nullptr;
23
24// See header for API docs.
25
26GHit::GHit(std::shared_ptr<GTouchable> gt,
27 const G4Step* thisStep,
28 const string& cScheme) :
29 G4VHit(),
30 colorSchema(cScheme),
31 gtouchable(gt) {
32 // Initialize per-step vectors if a step is provided.
33 if (thisStep) { addHitInfos(thisStep); }
34
35 // Cached derived quantities — computed lazily by calculateInfos() before publishing.
36 averageTime = UNINITIALIZEDNUMBERQUANTITY;
37 avgGlobalPosition = G4ThreeVector(UNINITIALIZEDNUMBERQUANTITY, UNINITIALIZEDNUMBERQUANTITY,
39 avgLocalPosition = G4ThreeVector(UNINITIALIZEDNUMBERQUANTITY, UNINITIALIZEDNUMBERQUANTITY,
41 processName = UNINITIALIZEDSTRINGQUANTITY;
42}
43
44bool GHit::is_same_hit(const GHit* hit) const {
45 if (!hit) // guard against nullptr
46 return false;
47
48 return *gtouchable == *(hit->getGTouchable());
49}
50
51vector<int> GHit::getTTID() const {
52 vector<int> ttid;
53 // Retrieve the identity vector from the associated GTouchable.
54 vector<GIdentifier> gids = getGID();
55 ttid.reserve(gids.size());
56 for (auto& gid : gids) {
57 // Push back the integer value of each identifier.
58 ttid.push_back(gid.getValue());
59 }
60 return ttid;
61}
62
63void GHit::Draw() {
64 auto visManager = G4VVisManager::GetConcreteInstance();
65 if (!visManager) return;
66
67 // Only care about schema if we are interactive.
68 setColorSchema();
69
70 // Check that globalPositions is not empty before accessing the first element.
71 if (globalPositions.empty()) return;
72
73 G4Circle circle(globalPositions[0]);
74 circle.SetFillStyle(G4Circle::filled);
75
76 double etot = getTotalEnergyDeposited();
77
78 if (etot > 0) {
79 circle.SetScreenSize(50);
80 circle.SetVisAttributes(G4VisAttributes(colour_hit));
81 }
82 else if (etot == 0) {
83 circle.SetScreenSize(15);
84 circle.SetVisAttributes(G4VisAttributes(colour_passby));
85 circle.SetFillStyle(G4Circle::hashed);
86 }
87
88 visManager->Draw(circle);
89}
90
91bool GHit::setColorSchema() {
92 // For now, hard-code the color schema.
93 colour_hit = G4Colour(1.0, 0.0, 0.0); // Red for hits with energy.
94 colour_passby = G4Colour(0.0, 1.0, 0.0); // Green for pass-by.
95 return false;
96}
97
99 trackVertexById.clear();
100 pdgById.clear();
101}
102
104 // This function is for testing purposes only.
105 // It randomizes the hit's global position and energy deposition.
106 // It should not be used in production code.
107
108 // Generate nsteps+1 entries to preserve the existing behavior exactly.
109 for (int i = 0; i < nsteps + 1; ++i) {
110 globalPositions.emplace_back(G4UniformRand() * 100, G4UniformRand() * 100, G4UniformRand() * 100);
111 localPositions.emplace_back(G4UniformRand() * 10, G4UniformRand() * 10, G4UniformRand() * 10);
112 trackVertexPositions.emplace_back(G4UniformRand() * 100, G4UniformRand() * 100, G4UniformRand() * 100);
113 motherTrackVertexPositions.emplace_back(UNINITIALIZEDNUMBERQUANTITY, UNINITIALIZEDNUMBERQUANTITY,
115 times.emplace_back(G4UniformRand() * 100);
116 edeps.emplace_back(G4UniformRand() * 10);
117 pids.emplace_back(11);
118 tids.emplace_back(i);
119 motherTids.emplace_back(0);
120 momenta.emplace_back(G4UniformRand() * 100, G4UniformRand() * 100, G4UniformRand() * 100);
121 trackEs.emplace_back(G4UniformRand() * 1000);
122 motherPids.emplace_back(2212); // proton as placeholder mother
123 processNames.emplace_back("placeholder");
124 }
125}
Stores step-by-step and aggregated information for a detector hit.
Definition ghit.h:38
static void clearTrackVertexCache()
Clear the per-thread track vertex cache.
Definition ghit.cc:98
size_t nsteps() const
Number of recorded steps.
Definition ghit.h:407
std::vector< GIdentifier > getGID() const
Get the detector element identity.
Definition ghit.h:483
std::vector< int > getTTID() const
Get the touchable identity values as integers.
Definition ghit.cc:51
bool is_same_hit(const GHit *hit) const
Compare this hit against another hit by sensitive-element identity.
Definition ghit.cc:44
void randomizeHitForTesting(int nsteps)
Randomize internal vectors for test-only usage.
Definition ghit.cc:103
GHit(std::shared_ptr< GTouchable > gt, const G4Step *thisStep=nullptr, const std::string &cScheme="default")
Construct a hit container and optionally seed it from a step.
Definition ghit.cc:26
void Draw() override
Visualize the hit using Geant4 visualization primitives.
Definition ghit.cc:63
double getTotalEnergyDeposited()
Get the total deposited energy across all recorded steps.
std::shared_ptr< GTouchable > getGTouchable() const
Get the associated sensitive-element descriptor.
Definition ghit.h:475
void addHitInfos(const G4Step *thisStep)
Append per-step information from a G4Step.
G4ThreadLocal G4Allocator< GHit > * GHitAllocator
Definition ghit.cc:22
#define UNINITIALIZEDNUMBERQUANTITY
#define UNINITIALIZEDSTRINGQUANTITY