19 #include "G4SDManager.hh"
20 #include "G4SystemOfUnits.hh"
21 #include "G4UnitsTable.hh"
32 :
GGSUserAction(), _messenger(this,
"/GGS/userActions/hitsAction/",
"Commands for hits action"), _outFileBase(
""),
33 _outFile(NULL), _outTree(NULL), _outTreeName(
""), _privateMessenger(this) {
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);
45 .SetGuidance(
"Class name of output integrated hits")
46 .SetStates(G4State_PreInit, G4State_Idle);
48 .SetGuidance(
"Class name of output particle hits")
49 .SetStates(G4State_PreInit, G4State_Idle);
51 .SetGuidance(
"Class name of output position hits")
52 .SetStates(G4State_PreInit, G4State_Idle);
60 for (
auto &handleElem : _outputHandles) {
61 OutputHandle &handle = handleElem.second;
62 if (handle.partHitArray) {
63 handle.partHitArray->SetOwner(
false);
65 if (handle.posHitArray) {
66 handle.posHitArray->SetOwner(
false);
74 static const std::string routineName(
"GSHitsAction::BeginOfRunAction");
80 if (!_outFile || (_outFile && _outFile->IsZombie())) {
81 COUT(ERROR) <<
"Cannot open ROOT file " << _outFileBase <<
"for run" << run->GetRunID() <<
ENDL;
88 if (_outTreeName ==
"") {
90 _outTree->SetTitle(TString(_outTree->GetTitle()) +
"Hits ");
92 _outTree =
new TTree(_outTreeName,
"GGS hits");
95 G4LogicalVolumeStore *logicalVS = G4LogicalVolumeStore::GetInstance();
99 std::vector<G4VSensitiveDetector *> processedSDs;
100 for (
unsigned int i = 0; i < logicalVS->size(); i++) {
101 G4LogicalVolume *logicalV = logicalVS->at(i);
104 std::vector<G4VSensitiveDetector *> sDets(0);
105 G4VSensitiveDetector *sd = logicalV->GetSensitiveDetector();
115 for (
auto &sdPtr : sDets) {
116 if (std::find(processedSDs.begin(), processedSDs.end(), sdPtr) != processedSDs.end()) {
119 auto sdName = sdPtr->GetName();
120 auto ggssdNamePos = sdName.find(
".GGSIntHitSD", 0);
121 if (ggssdNamePos != std::string::npos) {
123 if (intHitSD && intHitSD->isActive()) {
124 std::string sdNameReduced = sdName;
125 sdNameReduced.erase(ggssdNamePos, 12);
127 std::unordered_map<std::string, OutputHandle>::iterator handleIter = _outputHandles.find(sdName);
128 if (handleIter == _outputHandles.end()) {
131 auto insertResult = _outputHandles.insert(std::pair<std::string, OutputHandle>{sdName, OutputHandle()});
132 if (insertResult.second) {
133 handleIter = insertResult.first;
138 auto &newHandle = handleIter->second;
139 newHandle.intHitArray = std::unique_ptr<TClonesArray>(
new TClonesArray(
140 newHandle.intHitClassName.c_str()));
142 _outTree->Branch((sdNameReduced +
"-GGSIntHits").data(),
"TClonesArray", newHandle.intHitArray.get(), 64000);
146 auto thresholdIter = std::find_if(
147 _thresholds.begin(), _thresholds.end(), [&sdNameReduced](
const ThresholdStruct &threshold) {
148 return (threshold.hitType == 0) && (sdNameReduced == threshold.detectorName);
150 if (thresholdIter != _thresholds.end()) {
151 newHandle.intHitThreshold = thresholdIter->value;
155 newHandle.partHitArray = std::unique_ptr<TObjArray>(
new TObjArray());
158 newHandle.partHitArray->SetOwner(
true);
159 _outTree->Branch((sdNameReduced +
"-GGSPartHits").data(),
"TObjArray", newHandle.partHitArray.get(), 64000);
160 newHandle.partHitCAService =
163 auto thresholdIter = std::find_if(
164 _thresholds.begin(), _thresholds.end(), [&sdNameReduced](
const ThresholdStruct &threshold) {
165 return (threshold.hitType == 1) && (sdNameReduced == threshold.detectorName);
167 if (thresholdIter != _thresholds.end()) {
168 newHandle.partHitThreshold = thresholdIter->value;
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 =
178 auto thresholdIter = std::find_if(
179 _thresholds.begin(), _thresholds.end(), [&sdNameReduced](
const ThresholdStruct &threshold) {
180 return (threshold.hitType == 2) && (sdNameReduced == threshold.detectorName);
182 if (thresholdIter != _thresholds.end()) {
183 newHandle.posHitThreshold = thresholdIter->value;
188 unsigned int nBins = intHitSD->
GetTimeBins().size();
189 TArrayF timeBins(nBins);
190 for (
unsigned int iBin = 0; iBin < nBins; iBin++)
192 _outFile->WriteObject(&timeBins, (sdNameReduced +
"-GGSIntHits-TimeBins").data());
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());
203 processedSDs.push_back(sdPtr);
219 if (_outTreeName !=
"") {
237 _GiveBackAllArrays();
238 for (
auto &handleElem : _outputHandles) {
239 if (handleElem.second.posHitCAService) {
240 handleElem.second.posHitCAService->DeleteAll();
242 if (handleElem.second.partHitCAService) {
243 handleElem.second.partHitCAService->DeleteAll();
253 _GiveBackAllArrays();
255 for (
auto &handleElem : _outputHandles) {
256 OutputHandle &handle = handleElem.second;
257 if (handle.partHitArray) {
263 handle.partHitArray->SetOwner(
false);
265 handle.partHitArray->Clear();
268 handle.partHitArray->SetOwner(
true);
271 if (handle.posHitArray) {
272 handle.posHitArray->SetOwner(
false);
273 handle.posHitArray->Clear();
274 handle.posHitArray->SetOwner(
true);
288 G4SDManager *SDman = G4SDManager::GetSDMpointer();
291 G4HCofThisEvent *evHitColl =
event->GetHCofThisEvent();
294 for (
auto &handleElem : _outputHandles) {
296 const std::string &detectorName = handleElem.first;
297 OutputHandle &handle = handleElem.second;
299 G4int hitsCollID = SDman->GetCollectionID(detectorName);
302 handle.intHitArray->Clear();
303 _Convert(*((
GGSIntHitsCollection *)(evHitColl->GetHC(hitsCollID))), *(handle.intHitArray), handle, detectorName);
306 for (
int iIntHit = 0; iIntHit < handle.intHitArray->GetEntries(); ++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++) {
317 partHit->_posHitIndex = handle.posHitArray->GetEntries();
318 handle.posHitArray->AddAtAndExpand(partHit->_posHits, handle.posHitArray->GetEntries());
328 if (_outTreeName !=
"") {
336 void GGSHitsAction::_Convert(
const GGSPosHitsCollection &posHits, TClonesArray &tPosHits, OutputHandle &handle) {
338 if (tPosHits.GetEntries() != 0) {
339 throw std::runtime_error(
"GGSHitsAction::_Convert: dirty TClonesArray for position hits");
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) {
345 tPosHit.
eDep = posHit.GetEnergyDeposit() / GeV;
346 tPosHit.UserConversion(posHit);
355 void GGSHitsAction::_Convert(
const GGSPartHitsCollection &partHits, TClonesArray &tPartHits, OutputHandle &handle) {
357 if (tPartHits.GetEntries() != 0) {
358 throw std::runtime_error(
"GGSHitsAction::_Convert: dirty TClonesArray for particle hits");
360 for (
size_t iHit = 0, currTPartHit = 0; iHit < partHits.GetSize(); ++iHit) {
361 auto &partHit = *((
GGSPartHit *)(partHits.GetHit(iHit)));
363 auto &tPartHit = *((
GGSTPartHitBase *)(tPartHits.ConstructedAt(currTPartHit,
"C")));
366 if (!(tPartHit._posHits)) {
367 tPartHit._posHits = handle.posHitCAService->ProvideArray(
368 posHitColl->GetSize());
372 throw std::logic_error(
373 "GGSHitsAction::_Convert: particle hit already has position hits, this should not happen.");
375 _Convert(*posHitColl, *(tPartHit._posHits), handle);
378 if (partHit.GetEnergyDeposit() > handle.partHitThreshold || (posHitColl && tPartHit._posHits->GetEntries() != 0)) {
379 tPartHit.eDep = partHit.GetEnergyDeposit() / GeV;
380 tPartHit.UserConversion(partHit);
384 if (tPartHit._posHits) {
385 handle.posHitCAService->TakeBackArray(tPartHit._posHits);
387 tPartHits.RemoveLast();
394 void GGSHitsAction::_Convert(
const GGSIntHitsCollection &intHits, TClonesArray &tIntHits, OutputHandle &handle,
395 const std::string &detectorName) {
397 if (tIntHits.GetEntries() != 0) {
398 throw std::runtime_error(
"GGSHitsAction::_Convert: dirty TClonesArray for particle hits");
400 for (
size_t iHit = 0, currTIntHit = 0; iHit < intHits.GetSize(); ++iHit) {
401 auto &intHit = *((
GGSIntHit *)(intHits.GetHit(iHit)));
403 auto &tIntHit = *((
GGSTIntHitBase *)(tIntHits.ConstructedAt(currTIntHit,
"C")));
405 if (!tIntHit._partHits) {
406 tIntHit._partHits = handle.partHitCAService->ProvideArray(
407 partHitColl->GetSize());
411 throw std::logic_error(
412 "GGSHitsAction::_Convert: integrated hit already has particle hits, this should not happen.");
414 _Convert(*partHitColl, *(tIntHit._partHits), handle);
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;
426 tIntHit.UserConversion(intHit);
430 intHit.GetAbsTranslation(), intHit.GetID());
434 if (tIntHit._partHits) {
435 handle.partHitCAService->TakeBackArray(tIntHit._partHits);
439 tIntHits.RemoveLast();
447 const std::string routineName(
"[GGSHitsAction::SetOutputIntHitClass] ");
450 TClass *intHitClass = TClass::GetClass(className.c_str());
452 throw std::runtime_error(routineName +
"Dictionary for class " + className +
" not found");
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;
463 handleIter->second.intHitClassName = className;
469 const std::string routineName(
"[GGSHitsAction::SetOutputPartHitClass] ");
472 TClass *partHitClass = TClass::GetClass(className.c_str());
474 throw std::runtime_error(routineName +
"Dictionary for class " + className +
" not found");
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;
485 handleIter->second.partHitClassName = className;
491 const std::string routineName(
"[GGSHitsAction::SetOutputPosHitClass] ");
494 TClass *posHitClass = TClass::GetClass(className.c_str());
496 throw std::runtime_error(routineName +
"Dictionary for class " + className +
" not found");
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;
507 handleIter->second.posHitClassName = className;
513 const std::string &valueStr,
const std::string &unit) {
515 const std::string routineName(
"GGSHitsAction::SetHitThreshold");
518 if (hitType ==
"integrated") {
520 }
else if (hitType ==
"particle") {
522 }
else if (hitType ==
"position") {
525 throw std::runtime_error(routineName +
": Unknown hit type " + hitType +
" for detector " + detectorName);
530 value = std::stof(valueStr);
531 }
catch (std::invalid_argument &exc) {
532 COUT(ERROR) <<
"Invalid threshold for " << hitType <<
" hits of detector " << detectorName <<
ENDL;
534 }
catch (std::out_of_range &exc) {
535 COUT(ERROR) <<
"Invalid value for threshold of " << hitType <<
" hits of detector " << detectorName <<
ENDL;
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");
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();
550 if (unitIter == unitsList.end()) {
551 throw std::runtime_error(routineName +
": unit " + unit +
" not found");
554 auto thresholdIter = std::find_if(
555 _thresholds.begin(), _thresholds.end(), [&detectorName, hitTypeNo](
const ThresholdStruct &threshold) {
556 return (threshold.hitType == hitTypeNo) && (threshold.detectorName == detectorName);
558 if (thresholdIter == _thresholds.end()) {
559 _thresholds.push_back(ThresholdStruct{detectorName, hitTypeNo, -1.});
560 thresholdIter = _thresholds.end() - 1;
564 dynamic_cast<GGSIntHitSD *
>(G4SDManager::GetSDMpointer()->FindSensitiveDetector(detectorName +
".GGSIntHitSD"));
566 throw std::invalid_argument(routineName +
": detector " + detectorName +
" not found");
571 thresholdIter->value = value * (*unitIter)->GetValue();
572 sd->SetIntHitThreshold(thresholdIter->value);
575 thresholdIter->value = value * (*unitIter)->GetValue();
576 sd->SetPartHitThreshold(thresholdIter->value);
579 thresholdIter->value = value * (*unitIter)->GetValue();
580 sd->SetPosHitThreshold(thresholdIter->value);
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);
596 void GGSHitsAction::GGSHitsActionPrivateMessenger::SetNewValue(G4UIcommand *command, G4String newValue) {
597 const std::string routineName(
"GGSHitsAction::GGSHitsActionPrivateMessenger::SetNewValue");
598 if (command == _setThreshCmd.get()) {
600 if (tokens.size() != 4) {
601 throw std::runtime_error(routineName +
": setHitTreshold expects 4 parameters, but " +
602 std::to_string(tokens.size()) +
" were given.");
605 _hitsAction->SetHitThreshold(tokens[0], tokens[1], tokens[2], tokens[3]);
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++) {
617 if (intHit->_partHits) {
619 for (
int iPartHit = 0; iPartHit < nPartHits; iPartHit++) {
621 if (partHit->_posHits) {
622 handle.posHitCAService->TakeBackArray(partHit->_posHits);
625 handle.partHitCAService->TakeBackArray(intHit->_partHits);
G4THitsCollection< GGSPosHit > GGSPosHitsCollection
Alias for G4 template hits collection for GGSPosHit.
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.
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.
GGSHitsAction()
Constructor.
void CloseFileForThisRun(const path &baseName)
Closes the ROOT output file.
~GGSHitsAction()
Destructor.
#define COUT(level)
Smart log macro. It writes on stdout only if the specified verbosity level is lesser than the maximum...
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 of GGS Particle Hit.
Float_t eDep
Total deposited energy.
Definition of GGS Position Hit.
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.
A multiple sensitive detector.
void BeginOfEventAction(const G4Event *event)
Actions executed at beginning of event.
const std::vector< G4double > & GetTimeBins()
Time bins getter.
void BeginOfRunAction(const G4Run *run)
Actions executed at beginning of run.
bool GetPartHitsStorage()
Return current status of particle hits storage (i.e. persistence).
G4THitsCollection< GGSIntHit > GGSIntHitsCollection
Alias for G4 template hits collection for GGSIntHit.
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.
static GGSRootFileService & GetInstance()
Get reference to GGSRootFileService unique instance.
G4THitsCollection< GGSPartHit > GGSPartHitsCollection
Alias for G4 template hits collection for GGSPartHit.
Mother class for user actions in GGS.
Definition of GGS Integrated Hit.
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.
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.
Base class for storing G4 integrated hits.