GGS(GenericGEANT4Simulation)Software  2.6.0
 All Data Structures Namespaces Files Functions Variables Typedefs Macros
Public Member Functions
GGSIntHitSD Class Reference

Sensitive detector class for integrated hits. More...

#include <GGSIntHitSD.h>

Inheritance diagram for GGSIntHitSD:
Inheritance graph
[legend]
Collaboration diagram for GGSIntHitSD:
Collaboration graph
[legend]

Public Member Functions

 GGSIntHitSD (G4String name)
 Constructor. More...
 
 ~GGSIntHitSD ()
 Destructor.
 
void Initialize (G4HCofThisEvent *hitCollection)
 Initializes the sensitive detector. More...
 
G4bool ProcessHits (G4Step *aStep, G4TouchableHistory *ROHist)
 The hit processing method. More...
 
void EndOfEvent (G4HCofThisEvent *hitCollection)
 Sets the hits collection, as required by specifications of mother class. More...
 
void SetTimeBin (G4double time)
 Time bins setter. More...
 
const std::vector< G4double > & GetTimeBins ()
 Time bins getter. More...
 
void SetPartHitsStorage (bool store)
 Turn on or off the storage (i.e. persistency) of particle hits. More...
 
bool GetPartHitsStorage ()
 Return current status of particle hits storage (i.e. persistence). More...
 
void SetPosHitsStorage (bool store)
 Turn on or off the storage (i.e. persistence) of position hits. More...
 
bool GetPosHitsStorage ()
 Return current status of particle hits storage (i.e. persistence). More...
 

Detailed Description

Sensitive detector class for integrated hits.

This sensitive detector takes care of summing up all energy releases inside the sensitive volume it manages, so that a single hit is produced for each single volume with non-zero energy release for each event.

Definition at line 28 of file GGSIntHitSD.h.

Constructor & Destructor Documentation

GGSIntHitSD::GGSIntHitSD ( G4String  name)

Constructor.

Parameters
nameThe sensitive detector's name.

Definition at line 22 of file GGSIntHitSD.cpp.

22  :
23  G4VSensitiveDetector(name), _intHitCollection(NULL), _timeBins(0), _intHitClass("GGSIntHit"), _partHitClass(
24  "GGSPartHit"), _posHitClass("GGSPosHit"), _storePartHits(false), _storePosHits(false) {
25  collectionName.insert(name);
26  _messenger = new GGSIntHitSDMessenger(this);
27 
28  auto params = name.substr(name.find_first_of('.') + 1); // Remove volume name
29  auto secondPointPos = params.find_first_of('.');
30  if (secondPointPos != std::string::npos) {
31  // Parameters have been specified so let's handle them
32  params = params.substr(secondPointPos + 1); // Remove GGSIntHitSD
33  if (params.size() > 0) {
34  // Retrieve hit classes
35  std::string intHitClassName, partHitClassName, posHitClassName;
36  auto start = params.find_first_of('[');
37  auto end = params.find_first_of(']');
38  if (start < std::string::npos && end != std::string::npos) {
39  intHitClassName = params.substr(start + 1, end - start - 1);
40  }
41  else {
42  throw std::runtime_error(
43  std::string("[GGSIntHitSD::GGSIntHitSD] Missing definition of custom integrated hit class in \"") + params
44  + ("\""));
45  }
46 
47  start = params.find_first_of('[', end + 1);
48  end = params.find_first_of(']', start + 1);
49  if (start != std::string::npos && end != std::string::npos) {
50  partHitClassName = params.substr(start + 1, end - start - 1);
51  }
52  else {
53  throw std::runtime_error(
54  std::string("[GGSIntHitSD::GGSIntHitSD] Missing definition of custom particle hit class in \"") + params
55  + ("\""));
56  }
57 
58  start = params.find_first_of('[', end + 1);
59  end = params.find_first_of(']', start + 1);
60  if (start != std::string::npos && end != std::string::npos) {
61  posHitClassName = params.substr(start + 1, end - start - 1);
62  }
63  else {
64  throw std::runtime_error(
65  std::string("[GGSIntHitSD::GGSIntHitSD] Missing definition of custom position hit class in \"") + params
66  + ("\""));
67  }
68 
69  if (intHitClassName.size() > 0) {
71  if (std::find(buildersList.begin(), buildersList.end(), intHitClassName) == buildersList.end()) {
72  throw std::runtime_error(
73  std::string("[GGSIntHitSD::GGSIntHitSD] Integrated hit class \"") + intHitClassName + "\" not found.");
74  }
75  _intHitClass = intHitClassName;
76  }
77  if (partHitClassName.size() > 0) {
79  if (std::find(buildersList.begin(), buildersList.end(), partHitClassName) == buildersList.end()) {
80  throw std::runtime_error(
81  std::string("[GGSIntHitSD::GGSIntHitSD] Particle hit class \"") + partHitClassName + "\" not found.");
82  }
83  _partHitClass = partHitClassName;
84  }
85  if (posHitClassName.size() > 0) {
87  if (std::find(buildersList.begin(), buildersList.end(), posHitClassName) == buildersList.end()) {
88  throw std::runtime_error(
89  std::string("[GGSIntHitSD::GGSIntHitSD] Position hit class \"") + posHitClassName + "\" not found.");
90  }
91  _posHitClass = posHitClassName;
92  }
93  }
94  }
95 }
const std::vector< std::string > & GetListOfRegisteredBuilders()
Gets a vector containing the names of available registered builders.
Definition: GGSFactory.h:118
static GGSFactory & GetInstance()
Getter method for singleton pointer.
Definition: GGSFactory.h:75
The integrated hit SD messenger class.

Member Function Documentation

void GGSIntHitSD::EndOfEvent ( G4HCofThisEvent *  hitCollection)

Sets the hits collection, as required by specifications of mother class.

Parameters
hitCollectionThe hits collection.

Definition at line 253 of file GGSIntHitSD.cpp.

253  {
254 
255  G4int IHID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
256  hitCollection->AddHitsCollection(IHID, _intHitCollection);
257 
258 }
bool GGSIntHitSD::GetPartHitsStorage ( )
inline

Return current status of particle hits storage (i.e. persistence).

Returns
true if particle hits are going to be stored.

Definition at line 110 of file GGSIntHitSD.h.

110  {
111  return _storePartHits;
112  }
bool GGSIntHitSD::GetPosHitsStorage ( )
inline

Return current status of particle hits storage (i.e. persistence).

Returns
true if particle hits are going to be stored.

Definition at line 137 of file GGSIntHitSD.h.

137  {
138  return _storePosHits;
139  }
const std::vector<G4double>& GGSIntHitSD::GetTimeBins ( )
inline

Time bins getter.

First element is the end time of first bin (which always starts at t=0). Last element is the ending time of last bin but one, which corresponds to start time of last bin (last bin always ends at the end of event).

Returns
Vector of time bins.

Definition at line 83 of file GGSIntHitSD.h.

83  {
84  return _timeBins;
85  }
void GGSIntHitSD::Initialize ( G4HCofThisEvent *  hitCollection)

Initializes the sensitive detector.

Parameters
hitCollectionUnused (needed by the interface).

Definition at line 105 of file GGSIntHitSD.cpp.

105  {
106 
107  // The hits collection is added to G4HCofThisEvent in EndOfEvent.
108  // It is automatically deleted at the end of the event so it must be recreated every time.
109  // However, G4Allocator should provide fast memory allocation for this purpose.
110  _intHitCollection = new GGSIntHitsCollection(SensitiveDetectorName, collectionName[0]);
111 }
G4THitsCollection< GGSIntHit > GGSIntHitsCollection
Alias for G4 template hits collection for GGSIntHit.
Definition: GGSIntHit.h:272
G4bool GGSIntHitSD::ProcessHits ( G4Step *  aStep,
G4TouchableHistory *  ROHist 
)

The hit processing method.

For each physical sensitive volume, this method sums up the total energy release during this step to the total releases of previous steps.

Parameters
aStepThe current step
ROHistUnused (needed by the interface).
Returns
true if energy deposit of current step is greater than 0.

Definition at line 115 of file GGSIntHitSD.cpp.

115  {
116 
117  // The G4TouchableHistory argument of this method is non-null only if a readout geometry has been
118  // defined. Since this feature is not supported, it is currently always null.
119  // TODO: investigate the possibility to define a readout geometry
120 
121  if (aStep->GetTotalEnergyDeposit() == 0. && !_storePartHits)
122  return false;
123 
124  G4TouchableHandle touchable = aStep->GetPreStepPoint()->GetTouchableHandle();
125  G4ThreeVector volPos = touchable->GetTranslation(); // Absolute volume translation
126  G4VPhysicalVolume *physVol = touchable->GetVolume();
127 
128  // Build the absolute touchable ID
129  // Definition: the multiplicity of a physical volume is defined as the number of that volumes inside a given mother logical volume.
130  // Definition: given a history, the multiplicity at depth i is the multiplicity of the volume in history at depth i
131  // Definition: given a history, the copy number at depth i is the copy number of the volume in history at depth i
132  // Being N_i and M_i the copy number and multiplicity at depth i, the absolute ID of the touchable at history depth 0 is defined as:
133  // ID = N_0 + M_0*N_1 + M_0*M_1*N_2 + ...
134  // At each step, the product M_0*M_1*N_2*... is stored in the variable totMult.
135 
136  int historyDepth = touchable->GetHistoryDepth();
137  int totMult = 1;
138  int touchableID = 0;
139  // Map to store the multiplicity of each physical volume
140  static boost::unordered_map<G4VPhysicalVolume*, int> multMap;
141 
142  for (int iDepth = 0; iDepth < historyDepth; iDepth++) { // This loop will go up in history up to the daughter of the world volume (i.e. iDepth == historyDepth - 1)
143  touchableID += totMult * touchable->GetCopyNumber(iDepth);
144  // Compute total multiplicity for next step
145  if (iDepth != historyDepth - 1) { // No need to compute the total multiplicity on last step
146  int currVolMult = 0;
147  if (!(touchable->GetVolume(iDepth)->IsReplicated())) {
148  // See if the multiplicity for the current non-replicated volume has been computed before
149  boost::unordered_map<G4VPhysicalVolume*, int>::iterator storedMult = multMap.find(touchable->GetVolume(iDepth));
150  if (storedMult == multMap.end()) {
151  // Compute the multiplicity for non-replicated volumes placed inside the current mother logical
152  G4LogicalVolume *currLogVol = touchable->GetVolume(iDepth)->GetLogicalVolume();
153  G4LogicalVolume *motherLogVol = touchable->GetVolume(iDepth)->GetMotherLogical();
154  int nDaughters = motherLogVol->GetNoDaughters();
155  for (int iDaughter = 0; iDaughter < nDaughters; iDaughter++) {
156  if (motherLogVol->GetDaughter(iDaughter)->GetLogicalVolume() == currLogVol)
157  currVolMult++;
158  }
159  // Store the computed multiplicity
160  multMap[touchable->GetVolume(iDepth)] = currVolMult;
161  }
162  else {
163  // Use the stored multiplicity
164  currVolMult = storedMult->second;
165  }
166  }
167  else {
168  // Use the multiplicity of the replicated volume
169  currVolMult = touchable->GetVolume(iDepth)->GetMultiplicity();
170  }
171  totMult *= currVolMult;
172 
173  }
174  }
175 
176  /* *** Add step to integrated hit *** */
177  GGSIntHit *intHit = nullptr;
178  bool intFound = false;
179  // Search for existing integrated hit for this volume
180  for (int iHit = 0; iHit < _intHitCollection->entries(); iHit++) {
181  intHit = (*(_intHitCollection))[iHit];
182  if (physVol == intHit->GetVolume() && touchableID == intHit->GetID() && volPos == intHit->GetAbsTranslation()) {
183  intFound = true;
184  break;
185  }
186  }
187  if (!intFound) {
188  // If we reach this point, then this is a new hit
189  intHit = GGSIntHitFactory::GetInstance().CreateObject(_intHitClass, _timeBins).release();
190  if (!intHit) {
191  throw std::runtime_error(
192  std::string("[GGSIntHitSD::ProcessHits] Can't create integrated hit of class ").append(_intHitClass));
193  }
194  // Set global hit properties
195  intHit->SetAbsTranslation(volPos);
196  intHit->SetID(touchableID);
197  intHit->SetVolume(physVol);
198  intHit->SetPartHitsStorage(_storePartHits);
199  intHit->SetPosHitsStorage(_storePosHits);
200  // Insert the new integrated hit into the hit collection
201  _intHitCollection->insert(intHit);
202  }
203  // Add the current step to the hit
204  intHit->AddStep(*aStep);
205 
206  if (_storePartHits) {
207  /* *** Add step to particle hit *** */
208  GGSPartHit *partHit = nullptr;
209  auto *partHitCollection = intHit->GetPartHits();
210  bool partFound = false;
211  // Search for existing particle hit for this volume and particle
212  for (int iHit = 0; iHit < partHitCollection->entries(); iHit++) {
213  partHit = (*(partHitCollection))[iHit];
214  if ((*partHitCollection)[iHit]->GetTrackID() == aStep->GetTrack()->GetTrackID()) {
215  partFound = true;
216  break;
217  }
218  }
219  if (!partFound) {
220  // If we reach this point, then this is a new hit
221  partHit = GGSPartHitFactory::GetInstance().CreateObject(_partHitClass).release();
222  if (!partHit) {
223  throw std::runtime_error(
224  std::string("[GGSIntHitSD::ProcessHits] Can't create particle hit of class ").append(_partHitClass));
225  }
226 
227  // Set global hit properties
228  partHit->SetPosHitsStorage(_storePosHits);
229  // Insert the new particle hit into the hit collection
230  partHitCollection->insert(partHit);
231  }
232  // Add the current step to the hit
233  partHit->AddStep(*aStep);
234 
235  if (_storePosHits) {
236  /* *** Create position hit for current step *** */
237  GGSPosHit *posHit = GGSPosHitFactory::GetInstance().CreateObject(_posHitClass).release();
238  if (!posHit) {
239  throw std::runtime_error(
240  std::string("[GGSIntHitSD::ProcessHits] Can't create position hit of class ").append(_posHitClass));
241  }
242  // Insert the new position hit into the hit collection
243  partHit->GetPosHits()->insert(posHit);
244  // Set members of the position hit
245  posHit->SetStep(*aStep);
246  }
247  }
248 
249  return true;
250 }
virtual void AddStep(const G4Step &step)
Adds a step to the particle hit.
Definition: GGSIntHit.cpp:104
const G4ThreeVector & GetAbsTranslation()
Getter of absolute position of sensitive element.
Definition: GGSIntHit.h:242
virtual void AddStep(const G4Step &step)
Adds a step to the particle hit.
Definition: GGSPartHit.cpp:110
GGSPartHitsCollection * GetPartHits()
Getter of container of particle hits.
Definition: GGSIntHit.h:234
std::unique_ptr< T > CreateObject(const std::string &name, ConstructorArgs...arguments)
Create an object.
Definition: GGSFactory.h:123
void SetAbsTranslation(const G4ThreeVector &pos)
Set the position in global coordinates of the sensitive element.
Definition: GGSIntHit.h:136
void SetPosHitsStorage(bool flag)
Turn on or off the storage of position hits.
Definition: GGSIntHit.cpp:125
const G4VPhysicalVolume * GetVolume()
Getter for volume.
Definition: GGSIntHit.h:224
Definition of GGS Particle Hit.
Definition: GGSPartHit.h:35
Definition of GGS Position Hit.
Definition: GGSPosHit.h:31
virtual void SetStep(const G4Step &step)
Definition: GGSPosHit.cpp:65
G4int GetID()
Getter for volume ID.
Definition: GGSIntHit.h:250
void SetID(G4int id)
Setter for volume ID.
Definition: GGSIntHit.h:170
void SetPosHitsStorage(bool flag)
Turn on or off the storage of position hits.
Definition: GGSPartHit.cpp:130
void SetVolume(const G4VPhysicalVolume *volume)
Setter for volume.
Definition: GGSIntHit.h:99
GGSPosHitsCollection * GetPosHits()
Getter of container of position hits.
Definition: GGSPartHit.h:235
Definition of GGS Integrated Hit.
Definition: GGSIntHit.h:31
static GGSFactory & GetInstance()
Getter method for singleton pointer.
Definition: GGSFactory.h:75
void SetPartHitsStorage(bool flag)
Turn on or off the storage of particle hits.
Definition: GGSIntHit.cpp:144
void GGSIntHitSD::SetPartHitsStorage ( bool  store)
inline

Turn on or off the storage (i.e. persistency) of particle hits.

This method tells to the SD whether to activate the storage of particle hits or not. By default, storage is not activated.

Calling this method with "false" parameter will automatically disable position hits storage, i.e., it will call SetPartsHitStorage(false).

Parameters
storeif true, steps will be stored as position hits.
See Also
SetPosHitsStorage
GGSIntHit::SetPosHitsStorage

Definition at line 99 of file GGSIntHitSD.h.

99  {
100  _storePartHits = store;
101  if (store == false) {
102  SetPosHitsStorage(false);
103  }
104  }
void SetPosHitsStorage(bool store)
Turn on or off the storage (i.e. persistence) of position hits.
Definition: GGSIntHitSD.h:126
void GGSIntHitSD::SetPosHitsStorage ( bool  store)
inline

Turn on or off the storage (i.e. persistence) of position hits.

This method tells to the SD whether to activate the storage of position hits or not. By default, storage is not activated.

Calling this method with "true" parameter will automatically enable particle hits storage, i.e., it will call SetPartHitsStorage(true).

Parameters
storeif true, steps will be stored as position hits.
See Also
SetPartHitsStorage
GGSIntHit::SetPosHitsStorage

Definition at line 126 of file GGSIntHitSD.h.

126  {
127  _storePosHits = store;
128  if (store == true) {
129  SetPartHitsStorage(true);
130  }
131  }
void SetPartHitsStorage(bool store)
Turn on or off the storage (i.e. persistency) of particle hits.
Definition: GGSIntHitSD.h:99
void GGSIntHitSD::SetTimeBin ( G4double  time)
inline

Time bins setter.

This method adds a division in the timeline, effectively adding a time bin.

Parameters
timeThe time of first energy deposit.

Definition at line 70 of file GGSIntHitSD.h.

70  {
71  _timeBins.push_back(time);
72  std::sort(_timeBins.begin(), _timeBins.end());
73  }

The documentation for this class was generated from the following files: