19 #include "G4SystemOfUnits.hh"
20 #include "G4UnitsTable.hh"
21 #include "G4SDManager.hh"
32 GGSUserAction(), _messenger(this,
"/GGS/userActions/hitsAction/",
"Commands for hits action"), _outFileBase(
""), _outFile(
33 NULL), _outTree(NULL), _outTreeName(
""), _privateMessenger(this) {
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,
39 _messenger.DeclareProperty(
"treeName", _outTreeName).SetGuidance(
40 "Set the name of the TTree object in the ROOT output file.").SetStates(G4State_PreInit, G4State_Idle);
43 "Class name of output integrated hits").SetStates(G4State_PreInit, G4State_Idle);
45 "Class name of output particle hits").SetStates(G4State_PreInit, G4State_Idle);
47 "Class name of output position hits").SetStates(G4State_PreInit, G4State_Idle);
55 for (
auto &handleElem : _outputHandles) {
56 OutputHandle &handle = handleElem.second;
57 if (handle.partHitArray) {
58 handle.partHitArray->SetOwner(
false);
60 if (handle.posHitArray) {
61 handle.posHitArray->SetOwner(
false);
69 static const std::string routineName(
"GSHitsAction::BeginOfRunAction");
75 if (!_outFile || (_outFile && _outFile->IsZombie())) {
76 COUT(ERROR) <<
"Cannot open ROOT file " << _outFileBase <<
"for run" << run->GetRunID() <<
ENDL;
83 if (_outTreeName ==
"") {
85 _outTree->SetTitle(TString(_outTree->GetTitle()) +
"Hits ");
88 _outTree =
new TTree(_outTreeName,
"GGS hits");
91 G4LogicalVolumeStore *logicalVS = G4LogicalVolumeStore::GetInstance();
95 for (
unsigned int i = 0; i < logicalVS->size(); i++) {
96 G4LogicalVolume *logicalV = logicalVS->at(i);
99 std::vector<G4VSensitiveDetector*> sDets(0);
100 G4VSensitiveDetector *sd = logicalV->GetSensitiveDetector();
111 for (
auto &sdPtr : sDets) {
112 auto sdName = sdPtr->GetName();
113 auto ggssdNamePos = sdName.find(
".GGSIntHitSD", 0);
114 if (ggssdNamePos != std::string::npos) {
116 if (intHitSD && intHitSD->isActive()) {
117 std::string sdNameReduced = sdName;
118 sdNameReduced.erase(ggssdNamePos, 12);
120 std::unordered_map<std::string, OutputHandle>::iterator handleIter = _outputHandles.find(sdName);
121 if (handleIter == _outputHandles.end()) {
124 auto insertResult = _outputHandles.insert(std::pair<std::string, OutputHandle> { sdName, OutputHandle() });
125 if (insertResult.second) {
126 handleIter = insertResult.first;
131 auto &newHandle = handleIter->second;
132 newHandle.intHitArray = std::unique_ptr<TClonesArray>(
new TClonesArray(newHandle.intHitClassName.c_str()));
134 _outTree->Branch((sdNameReduced +
"-GGSIntHits").data(),
"TClonesArray", newHandle.intHitArray.get(), 64000);
138 auto thresholdIter = std::find_if(_thresholds.begin(), _thresholds.end(),
139 [&sdNameReduced](
const ThresholdStruct &threshold) {
140 return (threshold.hitType == 0) && (sdNameReduced == threshold.detectorName);
142 if (thresholdIter != _thresholds.end()) {
143 newHandle.intHitThreshold = thresholdIter->value;
147 newHandle.partHitArray = std::unique_ptr<TObjArray>(
new TObjArray());
150 newHandle.partHitArray->SetOwner(
true);
151 _outTree->Branch((sdNameReduced +
"-GGSPartHits").data(),
"TObjArray", newHandle.partHitArray.get(), 64000);
152 newHandle.partHitCAService = std::unique_ptr<GGSTClonesArrayService>(
155 auto thresholdIter = std::find_if(_thresholds.begin(), _thresholds.end(),
156 [&sdNameReduced](
const ThresholdStruct &threshold) {
157 return (threshold.hitType == 1) && (sdNameReduced == threshold.detectorName);
159 if (thresholdIter != _thresholds.end()) {
160 newHandle.partHitThreshold = thresholdIter->value;
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>(
170 auto thresholdIter = std::find_if(_thresholds.begin(), _thresholds.end(),
171 [&sdNameReduced](
const ThresholdStruct &threshold) {
172 return (threshold.hitType == 2) && (sdNameReduced == threshold.detectorName);
174 if (thresholdIter != _thresholds.end()) {
175 newHandle.posHitThreshold = thresholdIter->value;
180 unsigned int nBins = intHitSD->
GetTimeBins().size();
181 TArrayF timeBins(nBins);
182 for (
unsigned int iBin = 0; iBin < nBins; iBin++)
184 _outFile->WriteObject(&timeBins, (sdNameReduced +
"-GGSIntHits-TimeBins").data());
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());
211 if (_outTreeName !=
"") {
229 _GiveBackAllArrays();
230 for (
auto &handleElem : _outputHandles) {
231 if (handleElem.second.posHitCAService) {
232 handleElem.second.posHitCAService->DeleteAll();
234 if (handleElem.second.partHitCAService) {
235 handleElem.second.partHitCAService->DeleteAll();
245 _GiveBackAllArrays();
247 for (
auto &handleElem : _outputHandles) {
248 OutputHandle &handle = handleElem.second;
249 if (handle.partHitArray) {
255 handle.partHitArray->SetOwner(
false);
257 handle.partHitArray->Clear();
260 handle.partHitArray->SetOwner(
true);
263 if (handle.posHitArray) {
264 handle.posHitArray->SetOwner(
false);
265 handle.posHitArray->Clear();
266 handle.posHitArray->SetOwner(
true);
280 G4SDManager *SDman = G4SDManager::GetSDMpointer();
283 G4HCofThisEvent *evHitColl =
event->GetHCofThisEvent();
286 for (
auto &handleElem : _outputHandles) {
288 const std::string &detectorName = handleElem.first;
289 OutputHandle &handle = handleElem.second;
291 G4int hitsCollID = SDman->GetCollectionID(detectorName);
294 handle.intHitArray->Clear();
295 _Convert(*((
GGSIntHitsCollection*) (evHitColl->GetHC(hitsCollID))), *(handle.intHitArray), handle, detectorName);
298 for (
int iIntHit = 0; iIntHit < handle.intHitArray->GetEntries(); ++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++) {
309 partHit->_posHitIndex = handle.posHitArray->GetEntries();
310 handle.posHitArray->AddAtAndExpand(partHit->_posHits, handle.posHitArray->GetEntries());
321 if (_outTreeName !=
"") {
330 void GGSHitsAction::_Convert(
const GGSPosHitsCollection &posHits, TClonesArray &tPosHits, OutputHandle &handle) {
332 if (tPosHits.GetEntries() != 0) {
333 throw std::runtime_error(
"GGSHitsAction::_Convert: dirty TClonesArray for position hits");
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) {
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;
361 void GGSHitsAction::_Convert(
const GGSPartHitsCollection &partHits, TClonesArray &tPartHits, OutputHandle &handle) {
363 if (tPartHits.GetEntries() != 0) {
364 throw std::runtime_error(
"GGSHitsAction::_Convert: dirty TClonesArray for particle hits");
366 for (
size_t iHit = 0, currTPartHit = 0; iHit < partHits.GetSize(); ++iHit) {
367 auto &partHit = *((
GGSPartHit*) (partHits.GetHit(iHit)));
369 auto &tPartHit = *((
GGSTPartHit*) (tPartHits.ConstructedAt(currTPartHit,
"C")));
372 if (!(tPartHit._posHits)) {
373 tPartHit._posHits = handle.posHitCAService->ProvideArray(posHitColl->GetSize());
378 throw std::logic_error(
379 "GGSHitsAction::_Convert: particle hit already has position hits, this should not happen.");
381 _Convert(*posHitColl, *(tPartHit._posHits), handle);
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);
406 if (tPartHit._posHits) {
407 handle.posHitCAService->TakeBackArray(tPartHit._posHits);
410 tPartHits.RemoveLast();
418 void GGSHitsAction::_Convert(
const GGSIntHitsCollection &intHits, TClonesArray &tIntHits, OutputHandle &handle,
419 const std::string &detectorName) {
421 if (tIntHits.GetEntries() != 0) {
422 throw std::runtime_error(
"GGSHitsAction::_Convert: dirty TClonesArray for particle hits");
424 for (
size_t iHit = 0, currTIntHit = 0; iHit < intHits.GetSize(); ++iHit) {
425 auto &intHit = *((
GGSIntHit*) (intHits.GetHit(iHit)));
427 auto &tIntHit = *((
GGSTIntHit*) (tIntHits.ConstructedAt(currTIntHit,
"C")));
429 if (!tIntHit._partHits) {
430 tIntHit._partHits = handle.partHitCAService->ProvideArray(partHitColl->GetSize());
435 throw std::logic_error(
436 "GGSHitsAction::_Convert: integrated hit already has particle hits, this should not happen.");
438 _Convert(*partHitColl, *(tIntHit._partHits), handle);
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;
450 tIntHit.time = intHit.GetTime() / ns;
451 tIntHit.UserConversion(intHit);
455 intHit.GetAbsTranslation(), intHit.GetID());
460 if (tIntHit._partHits) {
461 handle.partHitCAService->TakeBackArray(tIntHit._partHits);
465 tIntHits.RemoveLast();
474 const std::string routineName(
"[GGSHitsAction::SetOutputIntHitClass] ");
477 TClass *intHitClass = TClass::GetClass(className.c_str());
479 throw std::runtime_error(routineName +
"Dictionary for class " + className +
" not found");
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;
490 handleIter->second.intHitClassName = className;
496 const std::string routineName(
"[GGSHitsAction::SetOutputPartHitClass] ");
499 TClass *partHitClass = TClass::GetClass(className.c_str());
501 throw std::runtime_error(routineName +
"Dictionary for class " + className +
" not found");
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;
512 handleIter->second.partHitClassName = className;
518 const std::string routineName(
"[GGSHitsAction::SetOutputPosHitClass] ");
521 TClass *posHitClass = TClass::GetClass(className.c_str());
523 throw std::runtime_error(routineName +
"Dictionary for class " + className +
" not found");
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;
534 handleIter->second.posHitClassName = className;
540 const std::string &valueStr,
const std::string &unit) {
542 const std::string routineName(
"GGSHitsAction::SetHitThreshold");
545 if (hitType ==
"integrated") {
548 else if (hitType ==
"particle") {
551 else if (hitType ==
"position") {
555 throw std::runtime_error(routineName +
": Unknown hit type " + hitType +
" for detector " + detectorName);
560 value = std::stof(valueStr);
561 }
catch (std::invalid_argument &exc) {
562 COUT(ERROR) <<
"Invalid threshold for " << hitType <<
" hits of detector " << detectorName <<
ENDL;
564 }
catch (std::out_of_range &exc) {
565 COUT(ERROR) <<
"Invalid value for threshold of " << hitType <<
" hits of detector " << detectorName <<
ENDL;
569 auto &unitsTable = G4UnitDefinition::GetUnitsTable();
570 auto unitCatIter = std::find_if(unitsTable.begin(), unitsTable.end(), [&unit](
const G4UnitsCategory *unitCat) {
571 return (unitCat->GetName() ==
"Energy");
573 if (unitCatIter == unitsTable.end()) {
574 throw std::runtime_error(routineName +
": energy unit category not found");
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();
581 if (unitIter == unitsList.end()) {
582 throw std::runtime_error(routineName +
": unit " + unit +
" not found");
585 auto thresholdIter = std::find_if(_thresholds.begin(), _thresholds.end(),
586 [&detectorName, hitTypeNo](
const ThresholdStruct &threshold) {
587 return (threshold.hitType == hitTypeNo) && (threshold.detectorName == detectorName);
589 if (thresholdIter == _thresholds.end()) {
590 _thresholds.push_back(ThresholdStruct { detectorName, hitTypeNo, -1. });
591 thresholdIter = _thresholds.end() - 1;
594 auto *sd =
dynamic_cast<GGSIntHitSD*
>(G4SDManager::GetSDMpointer()->FindSensitiveDetector(
595 detectorName +
".GGSIntHitSD"));
597 throw std::invalid_argument(routineName +
": detector " + detectorName +
" not found");
602 thresholdIter->value = value * (*unitIter)->GetValue();
603 sd->SetIntHitThreshold(thresholdIter->value);
606 thresholdIter->value = value * (*unitIter)->GetValue();
607 sd->SetPartHitThreshold(thresholdIter->value);
610 thresholdIter->value = value * (*unitIter)->GetValue();
611 sd->SetPosHitThreshold(thresholdIter->value);
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);
628 void GGSHitsAction::GGSHitsActionPrivateMessenger::SetNewValue(G4UIcommand *command, G4String newValue) {
629 const std::string routineName(
"GGSHitsAction::GGSHitsActionPrivateMessenger::SetNewValue");
630 if (command == _setThreshCmd.get()) {
632 if (tokens.size() != 4) {
633 throw std::runtime_error(
634 routineName +
": setHitTreshold expects 4 parameters, but " + std::to_string(tokens.size()) +
" were given.");
637 _hitsAction->SetHitThreshold(tokens[0], tokens[1], tokens[2], tokens[3]);
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++) {
650 if (intHit->_partHits) {
652 for (
int iPartHit = 0; iPartHit < nPartHits; iPartHit++) {
654 if (partHit->_posHits) {
655 handle.posHitCAService->TakeBackArray(partHit->_posHits);
658 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.
Float_t startMomentum[3]
Start momentum.
Float_t eDep
Total deposited energy.
Sensitive detector class for integrated hits.
Float_t time
Time of the hit.
virtual void UserConversion(const GGSPosHit &posHit)
Hook for user-defined conversion.
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.
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.
~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.
Class to store G4 particle hits.
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.
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.
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).
Int_t GetNPartHits()
Gets the number of particle hits.
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.
Class to store G4 position hits.
Class to store G4 position hits.
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.
Float_t startPoint[3]
Start point.
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.
Float_t endPoint[3]
End point.