ghit
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 
12 // MT definitions, as from:
13 // https://twiki.cern.ch/twiki/bin/view/Geant4/QuickMigrationGuideForGeant4V10
14 G4ThreadLocal G4Allocator<GHit>* GHitAllocator = nullptr;
15 
16 GHit::GHit(GTouchable *gt, const HitBitSet hbs, const G4Step *thisStep, string cScheme) : G4VHit(),
17  colorSchema(cScheme),
18  gtouchable(gt) {
19 
20  // initialize quantities based on HitBitSet, like globalPositions
21  if (thisStep) { addHitInfosForBitset(hbs, thisStep); }
22 
23  // unitialized quantities. To be calculated at the end of the steps by collectTrueInformation
24  // bit 0: always there
25  averageTime = UNINITIALIZEDNUMBERQUANTITY;
26  avgGlobalPosition = G4ThreeVector(UNINITIALIZEDNUMBERQUANTITY, UNINITIALIZEDNUMBERQUANTITY, UNINITIALIZEDNUMBERQUANTITY);
27  avgLocalPosition = G4ThreeVector(UNINITIALIZEDNUMBERQUANTITY, UNINITIALIZEDNUMBERQUANTITY, UNINITIALIZEDNUMBERQUANTITY);
28  processName = UNINITIALIZEDSTRINGQUANTITY;
29 
30 }
31 
32 bool GHit::is_same_hit(GHit *hit) {
33  return (*gtouchable) == (*hit->getGTouchable());
34 }
35 
36 
37 vector<int> GHit::getTTID() {
38  vector<int> ttid;
39  // Retrieve the identity vector from the associated GTouchable.
40  vector <GIdentifier> gids = getGID();
41  for (auto &gid: gids) {
42  // Push back the integer value of each identifier.
43  ttid.push_back(gid.getValue());
44  }
45  return ttid;
46 }
47 
48 
56  void GHit::Draw() {
57  G4VVisManager *pVVisManager = G4VVisManager::GetConcreteInstance();
58 
59  // only care about schema if we are interactive
60  if (pVVisManager) {
61  setColorSchema();
62 
63  // Check that globalPositions is not empty before accessing the first element.
64  if (globalPositions.empty()) return;
65 
66  G4Circle circle(globalPositions[0]);
67  circle.SetFillStyle(G4Circle::filled);
68 
69  double etot = getTotalEnergyDeposited();
70 
71  if (etot > 0) {
72  circle.SetScreenSize(10);
73  circle.SetVisAttributes(G4VisAttributes(colour_hit));
74  } else if (etot == 0) {
75  circle.SetScreenSize(8);
76  circle.SetVisAttributes(G4VisAttributes(colour_passby));
77  }
78 
79  pVVisManager->Draw(circle);
80  }
81 }
82 
83 
90 bool GHit::setColorSchema() {
91  // For now, hard-code the color schema.
92  colour_hit = G4Colour(1.0, 0.0, 0.0); // Red for hits with energy.
93  colour_passby = G4Colour(0.0, 1.0, 0.0); // Green for pass-by.
94  return false;
95 }
Represents a hit in the detector.
Definition: ghit.h:29
const vector< GIdentifier > getGID() const
Returns the detector element identity.
Definition: ghit.h:173
GHit(GTouchable *gt, const HitBitSet hbs, const G4Step *thisStep=nullptr, const string cScheme="default")
Constructor for GHit.
Definition: ghit.cc:16
void addHitInfosForBitset(const HitBitSet hbs, const G4Step *thisStep)
Adds hit information based on a HitBitSet.
Definition: addHitInfos.cc:17
bool is_same_hit(GHit *hit)
Compare two ghits.
Definition: ghit.cc:32
void Draw() override
Draws the hit using Geant4 visualization.
Definition: ghit.cc:56
float getTotalEnergyDeposited()
Computes the total energy deposited.
Definition: calculations.cc:76
vector< int > getTTID()
Returns the touchable identity values as integers.
Definition: ghit.cc:37
const GTouchable * getGTouchable() const
Returns the GTouchable associated with the hit.
Definition: ghit.h:167
std::bitset< NHITBITS > HitBitSet
G4ThreadLocal G4Allocator< GHit > * GHitAllocator
Definition: ghit.cc:14