GGS(GenericGEANT4Simulation)Software  2.6.3
 All Data Structures Namespaces Files Functions Variables Typedefs Macros
GGSHitsAction.cpp
Go to the documentation of this file.
1 /*
2  * GGSHitsAction.cpp
3  *
4  * Created on: Oct 8, 2010
5  * Author: Elena Vannuccini and Nicola Mori
6  */
7 
12 #include "utils/GGSNameDecoder.h"
13 #include "utils/GGSSmartLog.h"
17 #include "utils/GGSStringUtils.h"
18 
19 #include "G4SystemOfUnits.hh"
20 #include "G4UnitsTable.hh"
21 #include "G4SDManager.hh"
22 
23 #include "TArrayF.h"
24 
25 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
26 
28 
29 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
30 
32  GGSUserAction(), _messenger(this, "/GGS/userActions/hitsAction/", "Commands for hits action"), _outFileBase(""), _outFile(
33  NULL), _outTree(NULL), _outTreeName(""), _privateMessenger(this) {
34 
35  // Setup the messenger
36  _messenger.DeclareProperty("fileBase", _outFileBase).SetGuidance("Sets the base name for ROOT output file").SetGuidance(
37  " Can be with or without extension (.root will be used automatically)").SetStates(G4State_PreInit,
38  G4State_Idle);
39  _messenger.DeclareProperty("treeName", _outTreeName).SetGuidance(
40  "Set the name of the TTree object in the ROOT output file.").SetStates(G4State_PreInit, G4State_Idle);
41 
42  _messenger.DeclareMethod("outputIntHitClass", &GGSHitsAction::SetOutputIntHitClass).SetGuidance(
43  "Class name of output integrated hits").SetStates(G4State_PreInit, G4State_Idle);
44  _messenger.DeclareMethod("outputPartHitClass", &GGSHitsAction::SetOutputPartHitClass).SetGuidance(
45  "Class name of output particle hits").SetStates(G4State_PreInit, G4State_Idle);
46  _messenger.DeclareMethod("outputPosHitClass", &GGSHitsAction::SetOutputPosHitClass).SetGuidance(
47  "Class name of output position hits").SetStates(G4State_PreInit, G4State_Idle);
48 
49 }
50 
51 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
52 
54  // Let the TObjArry disown their TClonesArrays to avoid a double deletion since they are owned by the GGSTClonesArrayService.
55  for (auto &handleElem : _outputHandles) {
56  OutputHandle &handle = handleElem.second;
57  if (handle.partHitArray) {
58  handle.partHitArray->SetOwner(false);
59  }
60  if (handle.posHitArray) {
61  handle.posHitArray->SetOwner(false);
62  }
63  }
64 }
65 
66 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
67 
68 void GGSHitsAction::BeginOfRunAction(const G4Run *run) {
69  static const std::string routineName("GSHitsAction::BeginOfRunAction");
70 
71  // -----------------
72  // Open the run file
73  // -----------------
74  _outFile = GGSRootFileService::GetInstance().GetFileForThisRun(_outFileBase, run);
75  if (!_outFile || (_outFile && _outFile->IsZombie())) {
76  COUT(ERROR) << "Cannot open ROOT file " << _outFileBase << "for run" << run->GetRunID() << ENDL;
77  throw -1;
78  }
79 
80  // -----------------
81  // Create a new tree
82  // -----------------
83  if (_outTreeName == "") {
84  _outTree = GGSRootFileService::GetInstance().GetDefaultTree(_outFile);
85  _outTree->SetTitle(TString(_outTree->GetTitle()) + "Hits ");
86  }
87  else
88  _outTree = new TTree(_outTreeName, "GGS hits");
89 
90  // Loop over logical volumes...
91  G4LogicalVolumeStore *logicalVS = G4LogicalVolumeStore::GetInstance();
92  if (!logicalVS)
93  return;
94 
95  for (unsigned int i = 0; i < logicalVS->size(); i++) {
96  G4LogicalVolume *logicalV = logicalVS->at(i);
97 
98  // ...and select sensitive volumes with GGS hits SD
99  std::vector<G4VSensitiveDetector*> sDets(0);
100  G4VSensitiveDetector *sd = logicalV->GetSensitiveDetector();
101  if (sd) {
102  GGSMultiSensitiveDetector *msd = dynamic_cast<GGSMultiSensitiveDetector*>(sd);
103  if (msd) {
104  sDets = msd->GetListOfSensitiveDetectors();
105  }
106  else {
107  sDets.push_back(sd);
108  }
109  }
110 
111  for (auto &sdPtr : sDets) {
112  auto sdName = sdPtr->GetName();
113  auto ggssdNamePos = sdName.find(".GGSIntHitSD", 0);
114  if (ggssdNamePos != std::string::npos) {
115  GGSIntHitSD *intHitSD = dynamic_cast<GGSIntHitSD*>(sdPtr);
116  if (intHitSD && intHitSD->isActive()) {
117  std::string sdNameReduced = sdName;
118  sdNameReduced.erase(ggssdNamePos, 12); // Remove ".GGSIntHitSD"
119 
120  std::unordered_map<std::string, OutputHandle>::iterator handleIter = _outputHandles.find(sdName);
121  if (handleIter == _outputHandles.end()) {
122  // The handle for this detector does not exist yet, so create a standard one
123 
124  auto insertResult = _outputHandles.insert(std::pair<std::string, OutputHandle> { sdName, OutputHandle() });
125  if (insertResult.second) {
126  handleIter = insertResult.first;
127  }
128  }
129 
130  // Create arrays et al.
131  auto &newHandle = handleIter->second;
132  newHandle.intHitArray = std::unique_ptr<TClonesArray>(new TClonesArray(newHandle.intHitClassName.c_str())); // Create it here and Clear("C") it at beginning of event
133  // Create a branch to store hits in the output file
134  _outTree->Branch((sdNameReduced + "-GGSIntHits").data(), "TClonesArray", newHandle.intHitArray.get(), 64000);
135 
136  // Set output threshold
137  {
138  auto thresholdIter = std::find_if(_thresholds.begin(), _thresholds.end(),
139  [&sdNameReduced](const ThresholdStruct &threshold) {
140  return (threshold.hitType == 0) && (sdNameReduced == threshold.detectorName);
141  });
142  if (thresholdIter != _thresholds.end()) {
143  newHandle.intHitThreshold = thresholdIter->value;
144  }
145  }
146  if (intHitSD->GetPartHitsStorage()) {
147  newHandle.partHitArray = std::unique_ptr<TObjArray>(new TObjArray());
148  // Set ownership to avoid memory leaks when reading back the TObjArray.
149  // http://root.cern.ch/root/roottalk/roottalk02/0444.html
150  newHandle.partHitArray->SetOwner(true);
151  _outTree->Branch((sdNameReduced + "-GGSPartHits").data(), "TObjArray", newHandle.partHitArray.get(), 64000);
152  newHandle.partHitCAService = std::unique_ptr<GGSTClonesArrayService>(
153  new GGSTClonesArrayService(newHandle.partHitClassName.c_str()));
154  // Set output threshold
155  auto thresholdIter = std::find_if(_thresholds.begin(), _thresholds.end(),
156  [&sdNameReduced](const ThresholdStruct &threshold) {
157  return (threshold.hitType == 1) && (sdNameReduced == threshold.detectorName);
158  });
159  if (thresholdIter != _thresholds.end()) {
160  newHandle.partHitThreshold = thresholdIter->value;
161  }
162  }
163  if (intHitSD->GetPosHitsStorage()) {
164  newHandle.posHitArray = std::unique_ptr<TObjArray>(new TObjArray());
165  newHandle.posHitArray->SetOwner(true);
166  _outTree->Branch((sdNameReduced + "-GGSPosHits").data(), "TObjArray", newHandle.posHitArray.get(), 64000);
167  newHandle.posHitCAService = std::unique_ptr<GGSTClonesArrayService>(
168  new GGSTClonesArrayService(newHandle.posHitClassName.c_str()));
169  // Set output threshold
170  auto thresholdIter = std::find_if(_thresholds.begin(), _thresholds.end(),
171  [&sdNameReduced](const ThresholdStruct &threshold) {
172  return (threshold.hitType == 2) && (sdNameReduced == threshold.detectorName);
173  });
174  if (thresholdIter != _thresholds.end()) {
175  newHandle.posHitThreshold = thresholdIter->value;
176  }
177  }
178 
179  // Save the time bins as run products
180  unsigned int nBins = intHitSD->GetTimeBins().size();
181  TArrayF timeBins(nBins);
182  for (unsigned int iBin = 0; iBin < nBins; iBin++)
183  timeBins[iBin] = intHitSD->GetTimeBins()[iBin];
184  _outFile->WriteObject(&timeBins, (sdNameReduced + "-GGSIntHits-TimeBins").data());
185 
186  // Save the hit thresholds
187  std::vector<float> detThresholds { (float) (newHandle.intHitThreshold / GeV),
188  (float) (newHandle.partHitThreshold / GeV), (float) (newHandle.posHitThreshold / GeV) };
189  _outFile->WriteObject(&detThresholds, (sdNameReduced + "-GGSIntHits-Thresholds").data());
190 
191  } // Check on active volume
192 
193  } // Check on sensitive volume name
194 
195  } // End loop on sensitive detectors
196 
197  } // End loop on logical volumes
198 
199 }
200 
201 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
202 
203 void GGSHitsAction::EndOfRunAction(const G4Run*) {
204 
205 // -----------------
206 // close output file
207 // -----------------
208 
209  if (_outFile) {
210  _outFile->cd();
211  if (_outTreeName != "") {
212  // Write the tree if we are not using the default one
213  _outTree->Write();
214  }
215 
216  GGSRootFileService::GetInstance().CloseFileForThisRun(_outFileBase); // Will automatically save detector informations
217  _outFile = NULL;
218  _outTree = NULL; // The tree is delete by the file's destructor.
219  }
220  else {
221  throw -1; //TODO handle errors
222  }
223 
224 // -----------------
225 // delete detectors
226 // -----------------
227 
228  // Prepare for an eventual next run
229  _GiveBackAllArrays();
230  for (auto &handleElem : _outputHandles) {
231  if (handleElem.second.posHitCAService) {
232  handleElem.second.posHitCAService->DeleteAll();
233  }
234  if (handleElem.second.partHitCAService) {
235  handleElem.second.partHitCAService->DeleteAll();
236  }
237  }
238 }
239 
240 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
241 
242 void GGSHitsAction::BeginOfEventAction(const G4Event*) {
243 
244  // Give back arrays to the TCLonesArray services.
245  _GiveBackAllArrays();
246 
247  for (auto &handleElem : _outputHandles) {
248  OutputHandle &handle = handleElem.second;
249  if (handle.partHitArray) {
250 
251  // Clear the TObjarray persistence structures. These contains TClonesArrays which
252  // must not be deleted since they are owned and managed by the GGSTClonesArrayService...
253 
254  // ... so set the TObjArray to not own them...
255  handle.partHitArray->SetOwner(false);
256  // ... then clean to make size=0 without deleting the contained objects...
257  handle.partHitArray->Clear();
258  // ... and then restore ownership to avoid memory leak on readout of ROOT file.
259  // (see http://root.cern.ch/root/roottalk/roottalk02/0444.html)
260  handle.partHitArray->SetOwner(true);
261  }
262 
263  if (handle.posHitArray) {
264  handle.posHitArray->SetOwner(false);
265  handle.posHitArray->Clear();
266  handle.posHitArray->SetOwner(true);
267  }
268  }
269 }
270 
271 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
272 
273 void GGSHitsAction::EndOfEventAction(const G4Event *event) {
274 
275 //-------------------------------
276 // Set hits info
277 //-------------------------------
278 
279  // Retrieve SD manager
280  G4SDManager *SDman = G4SDManager::GetSDMpointer();
281 
282  // Retrieve event hits
283  G4HCofThisEvent *evHitColl = event->GetHCofThisEvent();
284 
285  // Loop on detectors
286  for (auto &handleElem : _outputHandles) {
287 
288  const std::string &detectorName = handleElem.first;
289  OutputHandle &handle = handleElem.second;
290  // Get hits of current detector
291  G4int hitsCollID = SDman->GetCollectionID(detectorName);
292 
293  // Convert hits to persistent data format
294  handle.intHitArray->Clear(); // Set entries to 0 but do not clear elements (will be cleared in _Convert)
295  _Convert(*((GGSIntHitsCollection*) (evHitColl->GetHC(hitsCollID))), *(handle.intHitArray), handle, detectorName);
296 
297  // Link particle and position hit collections to persistence array
298  for (int iIntHit = 0; iIntHit < handle.intHitArray->GetEntries(); ++iIntHit) {
299  GGSTIntHit *tIntHit = (GGSTIntHit*) (handle.intHitArray->At(iIntHit));
300  if (tIntHit->_partHits) {
301  handle.partHitArray->AddAtAndExpand(tIntHit->_partHits, handle.partHitArray->GetEntries());
302  if (tIntHit->_partHits->GetEntries() > 0 && ((GGSTPartHit*) (tIntHit->_partHits->At(0)))->_posHits) {
303  static TClonesArray *partHits = NULL;
304  static Int_t nPartHits = 0;
305  partHits = tIntHit->_partHits;
306  nPartHits = partHits->GetEntries();
307  for (Int_t iPartHit = 0; iPartHit < nPartHits; iPartHit++) {
308  GGSTPartHit *partHit = (GGSTPartHit*) (tIntHit->_partHits->At(iPartHit));
309  partHit->_posHitIndex = handle.posHitArray->GetEntries(); // Save reference to pos hit array
310  handle.posHitArray->AddAtAndExpand(partHit->_posHits, handle.posHitArray->GetEntries());
311  }
312  }
313  }
314 
315  }
316  }
317 
318 //-----------
319 // Fill tree
320 //-----------
321  if (_outTreeName != "") {
322  // Call fill only if not using the default tree
323  _outFile->cd();
324  _outTree->Fill();
325  }
326 
327 }
328 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
329 
330 void GGSHitsAction::_Convert(const GGSPosHitsCollection &posHits, TClonesArray &tPosHits, OutputHandle &handle) {
331 
332  if (tPosHits.GetEntries() != 0) {
333  throw std::runtime_error("GGSHitsAction::_Convert: dirty TClonesArray for position hits");
334  }
335  for (size_t iHit = 0, currTPosHit = 0; iHit < posHits.GetSize(); ++iHit) {
336  auto &posHit = *((GGSPosHit*) (posHits.GetHit(iHit)));
337  if (posHit.GetEnergyDeposit() > handle.posHitThreshold) {
338  GGSTPosHit &tPosHit = *((GGSTPosHit*) (tPosHits.ConstructedAt(currTPosHit, "C"))); // Get a cleared hit
339  tPosHit.eDep = posHit.GetEnergyDeposit() / GeV;
340  tPosHit.time = posHit.GetTime() / ns;
341  tPosHit.pathLength = posHit.GetPathLength() / cm;
342  tPosHit.startPoint[0] = posHit.GetStartPoint()[0] / cm;
343  tPosHit.startPoint[1] = posHit.GetStartPoint()[1] / cm;
344  tPosHit.startPoint[2] = posHit.GetStartPoint()[2] / cm;
345  tPosHit.endPoint[0] = posHit.GetEndPoint()[0] / cm;
346  tPosHit.endPoint[1] = posHit.GetEndPoint()[1] / cm;
347  tPosHit.endPoint[2] = posHit.GetEndPoint()[2] / cm;
348  tPosHit.startMomentum[0] = posHit.GetStartMomentum()[0] / GeV;
349  tPosHit.startMomentum[1] = posHit.GetStartMomentum()[1] / GeV;
350  tPosHit.startMomentum[2] = posHit.GetStartMomentum()[2] / GeV;
351  tPosHit.startEnergy = posHit.GetStartEnergy() / GeV;
352  tPosHit.UserConversion(posHit);
353  // Increase counter for output position hits
354  ++currTPosHit;
355  }
356  }
357 }
358 
359 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
360 
361 void GGSHitsAction::_Convert(const GGSPartHitsCollection &partHits, TClonesArray &tPartHits, OutputHandle &handle) {
362 
363  if (tPartHits.GetEntries() != 0) {
364  throw std::runtime_error("GGSHitsAction::_Convert: dirty TClonesArray for particle hits");
365  }
366  for (size_t iHit = 0, currTPartHit = 0; iHit < partHits.GetSize(); ++iHit) {
367  auto &partHit = *((GGSPartHit*) (partHits.GetHit(iHit)));
368  const GGSPosHitsCollection *posHitColl = partHit.GetPosHits();
369  auto &tPartHit = *((GGSTPartHit*) (tPartHits.ConstructedAt(currTPartHit, "C"))); // Get a cleared hit
370  if (posHitColl) {
371  // Convert the position hits
372  if (!(tPartHit._posHits)) {
373  tPartHit._posHits = handle.posHitCAService->ProvideArray(posHitColl->GetSize()); // Get a hit array with 0 entries and non-cleared elements
374  }
375  else {
376  // The particle hit should have given back its position hits array at BeginOfEvent, so this code
377  // should never be executed. Throw an exception.
378  throw std::logic_error(
379  "GGSHitsAction::_Convert: particle hit already has position hits, this should not happen.");
380  }
381  _Convert(*posHitColl, *(tPartHit._posHits), handle);
382  }
383  // Save the particle hit if its energy deposit is above threshold or if there are position hits above threshold
384  if (partHit.GetEnergyDeposit() > handle.partHitThreshold || (posHitColl && tPartHit._posHits->GetEntries() != 0)) {
385  tPartHit.eDep = partHit.GetEnergyDeposit() / GeV;
386  tPartHit.time = partHit.GetTime() / ns;
387  tPartHit.pathLength = partHit.GetPathLength() / cm;
388  tPartHit.entrancePoint[0] = partHit.GetEntrancePoint()[0] / cm;
389  tPartHit.entrancePoint[1] = partHit.GetEntrancePoint()[1] / cm;
390  tPartHit.entrancePoint[2] = partHit.GetEntrancePoint()[2] / cm;
391  tPartHit.exitPoint[0] = partHit.GetExitPoint()[0] / cm;
392  tPartHit.exitPoint[1] = partHit.GetExitPoint()[1] / cm;
393  tPartHit.exitPoint[2] = partHit.GetExitPoint()[2] / cm;
394  tPartHit.entranceMomentum[0] = partHit.GetEntranceMomentum()[0] / GeV;
395  tPartHit.entranceMomentum[1] = partHit.GetEntranceMomentum()[1] / GeV;
396  tPartHit.entranceMomentum[2] = partHit.GetEntranceMomentum()[2] / GeV;
397  tPartHit.entranceEnergy = partHit.GetEntranceEnergy() / GeV;
398  tPartHit.trackID = partHit.GetTrackID();
399  tPartHit.parentID = partHit.GetParentID();
400  tPartHit.particlePdg = partHit.GetParticlePdg();
401  tPartHit.UserConversion(partHit);
402  // Increase counter for output particle hits
403  ++currTPartHit;
404  }
405  else { // Otherwise remove the particle hit
406  if (tPartHit._posHits) { // Before deleting the particle hit we need to give back the array of position hits
407  handle.posHitCAService->TakeBackArray(tPartHit._posHits);
408 
409  }
410  tPartHits.RemoveLast();
411  }
412  }
413 
414 }
415 
416 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
417 
418 void GGSHitsAction::_Convert(const GGSIntHitsCollection &intHits, TClonesArray &tIntHits, OutputHandle &handle,
419  const std::string &detectorName) {
420 
421  if (tIntHits.GetEntries() != 0) {
422  throw std::runtime_error("GGSHitsAction::_Convert: dirty TClonesArray for particle hits");
423  }
424  for (size_t iHit = 0, currTIntHit = 0; iHit < intHits.GetSize(); ++iHit) {
425  auto &intHit = *((GGSIntHit*) (intHits.GetHit(iHit)));
426  const GGSPartHitsCollection *partHitColl = intHit.GetPartHits();
427  auto &tIntHit = *((GGSTIntHit*) (tIntHits.ConstructedAt(currTIntHit, "C"))); // Get a cleared hit
428  if (partHitColl) { // Convert the particle hits
429  if (!tIntHit._partHits) {
430  tIntHit._partHits = handle.partHitCAService->ProvideArray(partHitColl->GetSize()); // Get a hit array with 0 entries and non-cleared elements
431  }
432  else {
433  // The integrated hit should have given back its particle hits array at BeginOfEvent, so this code
434  // should never be executed. Throw an exception.
435  throw std::logic_error(
436  "GGSHitsAction::_Convert: integrated hit already has particle hits, this should not happen.");
437  }
438  _Convert(*partHitColl, *(tIntHit._partHits), handle);
439  }
440 
441  // Save the integrated hit if its energy deposit is above threshold or if there are particle hits
442  auto totEDep = intHit.GetTotalEnergyDep();
443  if (totEDep > handle.intHitThreshold || (partHitColl && tIntHit._partHits->GetEntries() != 0)) {
444  tIntHit.eDep = totEDep / GeV;
445  if (tIntHit.eDepTimeBin.GetSize() != (int) (intHit.GetEnergyDep().size()))
446  tIntHit.eDepTimeBin.Set(intHit.GetEnergyDep().size());
447  for (unsigned int i = 0; i < intHit.GetEnergyDep().size(); i++) {
448  tIntHit.eDepTimeBin[i] = intHit.GetEnergyDep()[i] / GeV;
449  }
450  tIntHit.time = intHit.GetTime() / ns;
451  tIntHit.UserConversion(intHit);
452 
453  // ----- Save volume informations of this hit ------
454  tIntHit._volumeIndex = GGSRootFileService::GetInstance().StoreVolume(_outFile, detectorName, intHit.GetVolume(),
455  intHit.GetAbsTranslation(), intHit.GetID());
456  // Increase counter for output integrated hits
457  ++currTIntHit;
458  }
459  else { // Otherwise remove the integrated hit
460  if (tIntHit._partHits) { // Before deleting the integrated hit we need to give back the array of particle hits
461  handle.partHitCAService->TakeBackArray(tIntHit._partHits);
462  // No need to give back arrays for position hits: this code portion is executed only when no particle hits are present,
463  // so there are no position hits
464  }
465  tIntHits.RemoveLast();
466  }
467  }
468 
469 }
470 
471 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
472 
473 void GGSHitsAction::SetOutputIntHitClass(std::string detectorName, const std::string &className) {
474  const std::string routineName("[GGSHitsAction::SetOutputIntHitClass] ");
475 
476  // Check if the given class has a dictionary
477  TClass *intHitClass = TClass::GetClass(className.c_str());
478  if (!intHitClass) {
479  throw std::runtime_error(routineName + "Dictionary for class " + className + " not found");
480  }
481 
482  // Check if the given detector already has a specified output class
483  detectorName.insert(detectorName.find_first_of('.'), ".GGSIntHitSD");
484  auto handleIter = _outputHandles.find(detectorName);
485  if (handleIter == _outputHandles.end()) {
486  handleIter = _outputHandles.insert(std::pair<std::string, OutputHandle> { detectorName, OutputHandle() }).first;
487  }
488 
489  // Set output facilities
490  handleIter->second.intHitClassName = className;
491 }
492 
493 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
494 
495 void GGSHitsAction::SetOutputPartHitClass(std::string detectorName, const std::string &className) {
496  const std::string routineName("[GGSHitsAction::SetOutputPartHitClass] ");
497 
498  // Check if the given class has a dictionary
499  TClass *partHitClass = TClass::GetClass(className.c_str());
500  if (!partHitClass) {
501  throw std::runtime_error(routineName + "Dictionary for class " + className + " not found");
502  }
503 
504  // Check if the given detector already has a specified output class
505  detectorName.insert(detectorName.find_first_of('.'), ".GGSIntHitSD");
506  auto handleIter = _outputHandles.find(detectorName);
507  if (handleIter == _outputHandles.end()) {
508  handleIter = _outputHandles.insert(std::pair<std::string, OutputHandle> { detectorName, OutputHandle() }).first;
509  }
510 
511  // Set output facilities
512  handleIter->second.partHitClassName = className;
513 }
514 
515 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
516 
517 void GGSHitsAction::SetOutputPosHitClass(std::string detectorName, const std::string &className) {
518  const std::string routineName("[GGSHitsAction::SetOutputPosHitClass] ");
519 
520  // Check if the given class has a dictionary
521  TClass *posHitClass = TClass::GetClass(className.c_str());
522  if (!posHitClass) {
523  throw std::runtime_error(routineName + "Dictionary for class " + className + " not found");
524  }
525 
526  // Check if the given detector already has a specified output class
527  detectorName.insert(detectorName.find_first_of('.'), ".GGSIntHitSD");
528  auto handleIter = _outputHandles.find(detectorName);
529  if (handleIter == _outputHandles.end()) {
530  handleIter = _outputHandles.insert(std::pair<std::string, OutputHandle> { detectorName, OutputHandle() }).first;
531  }
532 
533  // Set output facilities
534  handleIter->second.posHitClassName = className;
535 }
536 
537 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
538 
539 void GGSHitsAction::SetHitThreshold(const std::string &detectorName, const std::string &hitType,
540  const std::string &valueStr, const std::string &unit) {
541 
542  const std::string routineName("GGSHitsAction::SetHitThreshold");
543 
544  int hitTypeNo = -1;
545  if (hitType == "integrated") {
546  hitTypeNo = 0;
547  }
548  else if (hitType == "particle") {
549  hitTypeNo = 1;
550  }
551  else if (hitType == "position") {
552  hitTypeNo = 2;
553  }
554  else {
555  throw std::runtime_error(routineName + ": Unknown hit type " + hitType + " for detector " + detectorName);
556  }
557 
558  float value = 0.;
559  try {
560  value = std::stof(valueStr);
561  } catch (std::invalid_argument &exc) {
562  COUT(ERROR) << "Invalid threshold for " << hitType << " hits of detector " << detectorName << ENDL;
563  throw exc;
564  } catch (std::out_of_range &exc) {
565  COUT(ERROR) << "Invalid value for threshold of " << hitType << " hits of detector " << detectorName << ENDL;
566  throw exc;
567  }
568 
569  auto &unitsTable = G4UnitDefinition::GetUnitsTable();
570  auto unitCatIter = std::find_if(unitsTable.begin(), unitsTable.end(), [&unit](const G4UnitsCategory *unitCat) {
571  return (unitCat->GetName() == "Energy");
572  });
573  if (unitCatIter == unitsTable.end()) {
574  throw std::runtime_error(routineName + ": energy unit category not found");
575  }
576 
577  auto &unitsList = (*unitCatIter)->GetUnitsList();
578  auto unitIter = std::find_if(unitsList.begin(), unitsList.end(), [&unit](const G4UnitDefinition *unitDef) {
579  return unit == unitDef->GetSymbol().c_str();
580  });
581  if (unitIter == unitsList.end()) {
582  throw std::runtime_error(routineName + ": unit " + unit + " not found");
583  }
584 
585  auto thresholdIter = std::find_if(_thresholds.begin(), _thresholds.end(),
586  [&detectorName, hitTypeNo](const ThresholdStruct &threshold) {
587  return (threshold.hitType == hitTypeNo) && (threshold.detectorName == detectorName);
588  });
589  if (thresholdIter == _thresholds.end()) {
590  _thresholds.push_back(ThresholdStruct { detectorName, hitTypeNo, -1. });
591  thresholdIter = _thresholds.end() - 1;
592  }
593 
594  auto *sd = dynamic_cast<GGSIntHitSD*>(G4SDManager::GetSDMpointer()->FindSensitiveDetector(
595  detectorName + ".GGSIntHitSD"));
596  if (!sd) {
597  throw std::invalid_argument(routineName + ": detector " + detectorName + " not found");
598  }
599 
600  switch (hitTypeNo) {
601  case 0:
602  thresholdIter->value = value * (*unitIter)->GetValue();
603  sd->SetIntHitThreshold(thresholdIter->value);
604  break;
605  case 1:
606  thresholdIter->value = value * (*unitIter)->GetValue();
607  sd->SetPartHitThreshold(thresholdIter->value);
608  break;
609  case 2:
610  thresholdIter->value = value * (*unitIter)->GetValue();
611  sd->SetPosHitThreshold(thresholdIter->value);
612  break;
613  }
614 
615 }
616 
617 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
618 
619 GGSHitsAction::GGSHitsActionPrivateMessenger::GGSHitsActionPrivateMessenger(GGSHitsAction *hitsAction) :
620  _hitsAction(hitsAction) {
621  _setThreshCmd = std::make_unique<G4UIcmdWithAString>("/GGS/userActions/hitsAction/setHitThreshold", this);
622  _setThreshCmd->SetGuidance("Set energy deposit threshold for hit output");
623  _setThreshCmd->AvailableForStates(G4State_PreInit, G4State_Idle);
624 }
625 
626 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
627 
628 void GGSHitsAction::GGSHitsActionPrivateMessenger::SetNewValue(G4UIcommand *command, G4String newValue) {
629  const std::string routineName("GGSHitsAction::GGSHitsActionPrivateMessenger::SetNewValue");
630  if (command == _setThreshCmd.get()) {
631  auto tokens = GGSStringUtils::Tokenize(newValue);
632  if (tokens.size() != 4) {
633  throw std::runtime_error(
634  routineName + ": setHitTreshold expects 4 parameters, but " + std::to_string(tokens.size()) + " were given.");
635  }
636 
637  _hitsAction->SetHitThreshold(tokens[0], tokens[1], tokens[2], tokens[3]);
638 
639  }
640 }
641 
642 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
643 
644 void GGSHitsAction::_GiveBackAllArrays() {
645  for (auto &handleElem : _outputHandles) {
646  OutputHandle &handle = handleElem.second;
647  int nIntHits = handle.intHitArray->GetEntries();
648  for (int iIntHit = 0; iIntHit < nIntHits; iIntHit++) {
649  GGSTIntHit *intHit = (GGSTIntHit*) (handle.intHitArray->At(iIntHit));
650  if (intHit->_partHits) {
651  int nPartHits = intHit->GetNPartHits();
652  for (int iPartHit = 0; iPartHit < nPartHits; iPartHit++) {
653  GGSTPartHit *partHit = intHit->GetPartHit(iPartHit);
654  if (partHit->_posHits) {
655  handle.posHitCAService->TakeBackArray(partHit->_posHits);
656  }
657  }
658  handle.partHitCAService->TakeBackArray(intHit->_partHits);
659  }
660  }
661  }
662 }
G4THitsCollection< GGSPosHit > GGSPosHitsCollection
Alias for G4 template hits collection for GGSPosHit.
Definition: GGSPosHit.h:223
TFile * GetFileForThisRun(const path &baseName, const G4Run *run)
Opens a file for a given run and returns a pointer to it.
Float_t startMomentum[3]
Start momentum.
Definition: GGSTHits.h:42
Float_t eDep
Total deposited energy.
Definition: GGSTHits.h:37
Sensitive detector class for integrated hits.
Definition: GGSIntHitSD.h:28
Float_t time
Time of the hit.
Definition: GGSTHits.h:38
virtual void UserConversion(const GGSPosHit &posHit)
Hook for user-defined conversion.
Definition: GGSTHits.h:67
A service which manages a pool of reusable TClonesArray.
void EndOfRunAction(const G4Run *run)
Actions executed at end of run.
const std::vector< G4VSensitiveDetector * > & GetListOfSensitiveDetectors()
Retrieves a list of sensitive detectors collected into this object.
Float_t pathLength
Path length of the hit.
Definition: GGSTHits.h:39
void SetHitThreshold(const std::string &detectorName, const std::string &hitType, const std::string &valueStr, const std::string &unit)
Sets a lower energy deposit threshold for hits.
GGSHitsAction()
Constructor.
void CloseFileForThisRun(const path &baseName)
Closes the ROOT output file.
Float_t startEnergy
Start kinetic energy.
Definition: GGSTHits.h:43
#define ENDL
Definition: GGSSmartLog.h:104
~GGSHitsAction()
Destructor.
#define COUT(level)
Smart log macro. It writes on stdout only if the specified verbosity level is lesser than the maximum...
Definition: GGSSmartLog.h:77
void EndOfEventAction(const G4Event *event)
Actions executed at end of event.
Class to store G4 particle hits.
Definition: GGSTHits.h:82
Tokens Tokenize(const std::string &str, char delimiter= ' ')
Extracts words from a string.
bool GetPosHitsStorage()
Return current status of particle hits storage (i.e. persistence).
Definition: GGSIntHitSD.h:137
Definition of GGS Particle Hit.
Definition: GGSPartHit.h:35
Definition of GGS Position Hit.
Definition: GGSPosHit.h:31
int StoreVolume(TFile *filePtr, const std::string &detector, const G4VPhysicalVolume *volume, const G4ThreeVector &position, int id)
Set persistence for the specified volume.
A multiple sensitive detector.
void BeginOfEventAction(const G4Event *event)
Actions executed at beginning of event.
const std::vector< G4double > & GetTimeBins()
Time bins getter.
Definition: GGSIntHitSD.h:83
void BeginOfRunAction(const G4Run *run)
Actions executed at beginning of run.
bool GetPartHitsStorage()
Return current status of particle hits storage (i.e. persistence).
Definition: GGSIntHitSD.h:110
Int_t GetNPartHits()
Gets the number of particle hits.
Definition: GGSTHits.cpp:240
G4THitsCollection< GGSIntHit > GGSIntHitsCollection
Alias for G4 template hits collection for GGSIntHit.
Definition: GGSIntHit.h:272
void SetOutputPartHitClass(std::string detectorName, const std::string &className)
Sets the name of the particle hit class to be used for output.
TTree * GetDefaultTree(TFile *file)
Gets the default tree for this file.
Class to store G4 position hits.
Definition: GGSTHits.h:34
Class to store G4 position hits.
Definition: GGSTHits.h:160
An action which saves hits in sensitive volumes on ROOT file.
Definition: GGSHitsAction.h:63
static GGSRootFileService & GetInstance()
Get reference to GGSRootFileService unique instance.
G4THitsCollection< GGSPartHit > GGSPartHitsCollection
Alias for G4 template hits collection for GGSPartHit.
Definition: GGSPartHit.h:259
Mother class for user actions in GGS.
Definition: GGSUserAction.h:27
Definition of GGS Integrated Hit.
Definition: GGSIntHit.h:31
Float_t startPoint[3]
Start point.
Definition: GGSTHits.h:40
void SetOutputPosHitClass(std::string detectorName, const std::string &className)
Sets the name of the position hit class to be used for output.
#define RegisterUA(uaClassName)
Macro for registration of user actions classes.
void SetOutputIntHitClass(std::string detectorName, const std::string &className)
Sets the name of the integrated hit class to be used for output of a given detector.
GGSTPartHit * GetPartHit(unsigned int iHit)
Get the specified particle hit.
Definition: GGSTHits.cpp:249
Float_t endPoint[3]
End point.
Definition: GGSTHits.h:41