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