ghit
Loading...
Searching...
No Matches
ghit.cc
Go to the documentation of this file.
1// ghit
2#include "ghit.h"
3
4#include <utility>
5
6// glibrary
7#include "gutsConventions.h"
8
9// geant4
10#include "G4VVisManager.hh"
11#include "G4Circle.hh"
12#include "G4VisAttributes.hh"
13#include "Randomize.hh"
14
15using std::string;
16using std::vector;
17
18std::atomic<int> GHit::globalHitCounter{0};
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 HitBitSet hbs,
28 const G4Step* thisStep,
29 const string& cScheme) :
30 G4VHit(),
31 colorSchema(cScheme),
32 gtouchable(gt) {
33 // Initialize per-step vectors if a step is provided.
34 if (thisStep) { addHitInfosForBitset(hbs, thisStep); }
35
36 // Uninitialized quantities, to be calculated at the end of the steps by collectTrueInformation
37 // bit 0: always there
38 averageTime = UNINITIALIZEDNUMBERQUANTITY;
39 avgGlobalPosition = G4ThreeVector(UNINITIALIZEDNUMBERQUANTITY, UNINITIALIZEDNUMBERQUANTITY,
41 avgLocalPosition = G4ThreeVector(UNINITIALIZEDNUMBERQUANTITY, UNINITIALIZEDNUMBERQUANTITY,
43 processName = UNINITIALIZEDSTRINGQUANTITY;
44}
45
46bool GHit::is_same_hit(const GHit* hit) const {
47 if (!hit) // guard against nullptr
48 return false;
49
50 return *gtouchable == *(hit->getGTouchable());
51}
52
53vector<int> GHit::getTTID() const {
54 vector<int> ttid;
55 // Retrieve the identity vector from the associated GTouchable.
56 vector<GIdentifier> gids = getGID();
57 ttid.reserve(gids.size());
58 for (auto& gid : gids) {
59 // Push back the integer value of each identifier.
60 ttid.push_back(gid.getValue());
61 }
62 return ttid;
63}
64
65void GHit::Draw() {
66 auto visManager = G4VVisManager::GetConcreteInstance();
67 if (!visManager) return;
68
69 // Only care about schema if we are interactive.
70 setColorSchema();
71
72 // Check that globalPositions is not empty before accessing the first element.
73 if (globalPositions.empty()) return;
74
75 G4Circle circle(globalPositions[0]);
76 circle.SetFillStyle(G4Circle::filled);
77
78 double etot = getTotalEnergyDeposited();
79
80 if (etot > 0) {
81 circle.SetScreenSize(50);
82 circle.SetVisAttributes(G4VisAttributes(colour_hit));
83 }
84 else if (etot == 0) {
85 circle.SetScreenSize(15);
86 circle.SetVisAttributes(G4VisAttributes(colour_passby));
87 circle.SetFillStyle(G4Circle::hashed);
88 }
89
90 visManager->Draw(circle);
91}
92
93bool GHit::setColorSchema() {
94 // For now, hard-code the color schema.
95 colour_hit = G4Colour(1.0, 0.0, 0.0); // Red for hits with energy.
96 colour_passby = G4Colour(0.0, 1.0, 0.0); // Green for pass-by.
97 return false;
98}
99
101 // This function is for testing purposes only.
102 // It randomizes the hit's global position and energy deposition.
103 // It should not be used in production code.
104
105 // Generate nsteps+1 entries to preserve the existing behavior exactly.
106 for (int i = 0; i < nsteps + 1; ++i) {
107 globalPositions.emplace_back(G4UniformRand() * 100, G4UniformRand() * 100, G4UniformRand() * 100);
108 localPositions.emplace_back(G4UniformRand() * 10, G4UniformRand() * 10, G4UniformRand() * 10);
109 times.emplace_back(G4UniformRand() * 100);
110 edeps.emplace_back(G4UniformRand() * 10);
111 Es.emplace_back(G4UniformRand() * 10);
112
113 pids.emplace_back(static_cast<int>(G4UniformRand() * 1000)); // Random particle ID
114 }
115}
Stores step-by-step and aggregated information for a detector hit.
Definition ghit.h:43
GHit(std::shared_ptr< GTouchable > gt, HitBitSet hbs, 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 addHitInfosForBitset(HitBitSet hbs, const G4Step *thisStep)
Append per-step information from a G4Step according to a bitset.
size_t nsteps() const
Number of recorded steps for the optional-energy vector.
Definition ghit.h:322
std::vector< GIdentifier > getGID() const
Get the detector element identity.
Definition ghit.h:345
std::vector< int > getTTID() const
Get the touchable identity values as integers.
Definition ghit.cc:53
bool is_same_hit(const GHit *hit) const
Compare this hit against another hit by sensitive-element identity.
Definition ghit.cc:46
void randomizeHitForTesting(int nsteps)
Randomize internal vectors for test-only usage.
Definition ghit.cc:100
void Draw() override
Visualize the hit using Geant4 visualization primitives.
Definition ghit.cc:65
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:337
std::bitset< NHITBITS > HitBitSet
Bitset selecting which optional hit information is recorded.
G4ThreadLocal G4Allocator< GHit > * GHitAllocator
Definition ghit.cc:22
#define UNINITIALIZEDNUMBERQUANTITY
#define UNINITIALIZEDSTRINGQUANTITY