ATLAS Offline Software
Loading...
Searching...
No Matches
sTGCSensitiveDetector.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2025 CERN for the benefit of the ATLAS collaboration
3*/
4
8#include "G4Geantino.hh"
9#include "G4ChargedGeantino.hh"
10
11#include "G4Exception.hh"
12#include "G4Track.hh"
13
16
17#include <string>
18
19
20// construction/destruction
22 const std::string& hitCollectionName,
23 unsigned baseDepth):
24 G4VSensitiveDetector( name ),
25 AthMessaging{name},
26 m_hitCollectionName( hitCollectionName ),
27 m_baseDepth{baseDepth} {}
28
29// Implemenation of memebr functions
31{
32 m_sTGCSimHitCollection = nullptr;
33 if (auto* eventInfo = AtlasG4EventUserInfo::GetEventUserInfo()) {
34 m_sTGCSimHitCollection = eventInfo->GetHitCollectionMap()->Find<sTGCSimHitCollection>(m_hitCollectionName);
35 m_g4UserEventInfo = eventInfo;
36 }
37}
38
39G4bool sTGCSensitiveDetector::ProcessHits(G4Step* aStep,G4TouchableHistory* /*ROHist*/)
40{
42 G4Exception("sTGCSensitiveDetector::ProcessHits", "sTGCHitCollectionMissing", FatalException,
43 "Hit collection not initialized; did SetupEvent run?");
44 return false;
45 }
46 G4Track* currentTrack = aStep->GetTrack();
47 int charge=currentTrack->GetDefinition()->GetPDGCharge();
48
49 bool geantinoHit = (currentTrack->GetDefinition()==G4Geantino::GeantinoDefinition()) ||
50 (currentTrack->GetDefinition()==G4ChargedGeantino::ChargedGeantinoDefinition());
51
52 if (!charge && (!geantinoHit)) return false;
53 // G4cout << "\t\t sTGCSD: Hit in a sensitive layer!!!!! " << G4endl;
54 G4StepPoint* postStep=aStep->GetPostStepPoint();
55 const G4Step* post_Step=aStep->GetTrack()->GetStep();
56
57 Amg::Vector3D position = Amg::Hep3VectorToEigen(postStep->GetPosition());
58
59 G4StepPoint* preStep = aStep->GetPreStepPoint();
60 Amg::Vector3D preposition = Amg::Hep3VectorToEigen(preStep->GetPosition());
61
62 int pdgCode=currentTrack->GetDefinition()->GetPDGEncoding();
63
64 float globalTime=postStep->GetGlobalTime();
65
66 Amg::Vector3D direction= Amg::Hep3VectorToEigen( postStep->GetMomentumDirection() );
67 float depositEnergy=post_Step->GetTotalEnergyDeposit();
68
69 if (depositEnergy<0.0001 && (!geantinoHit)) return false;
70
71 const G4TouchableHistory* touchHist = static_cast<const G4TouchableHistory*>(aStep->GetPreStepPoint()->GetTouchable());
72
73 // int iDepth=touchHist->GetHistoryDepth();
74 // G4cout << "\t\t\t\t Touchable history dump "<<G4endl;
75 int nLayer=touchHist->GetVolume(m_baseDepth)->GetCopyNo();
76 std::string chName=touchHist->GetVolume(1+m_baseDepth)->GetLogicalVolume()->GetName();
77 //G4cout << "sTGCSensitiveDetector name: "<<chName<<G4endl;
78 std::string subType=chName.substr(chName.find('-')+1);
79 //G4cout << "\t\t sType: "<<subType);
80 if (subType[0]!='T'&&subType[0]!='Q' ) ATH_MSG_WARNING(" something is wrong, this is no sTGC! "<<chName<<", "<<Amg::toString(preposition));
81 int iRing = std::atoi(&subType[2]) -1;
82 ATH_MSG_VERBOSE("Volume name: "<<chName<<", nLayer: "<<nLayer<<", subType: "<<subType<<", iRing: "<<iRing);
83
84 // identifiers have eta naming 0-2, eta encoded in subtype is 1-3
85 // double phiDiff=2*M_PI;
86
87 G4ThreeVector posH=postStep->GetPosition(); //posH is equivalent to position - eigen not used to avoid additional dependence on EventPrimitives
88 if (subType[1]=='L') posH.rotateZ(M_PI/8.);
89 double phiHit=posH.phi();
90 if(phiHit<=0) phiHit+=2.*M_PI;
91 int iPhi=1+int(phiHit/(M_PI/4.));
92 iPhi*=2;
93 if (subType[1]=='L') iPhi-=1;
94
95 int iSide=1;
96 if (position.z()<0) iSide=-1;
97
98 int mLayer=0;
99 if (subType[1]=='S')
100 {
101 if (subType[3]=='C') mLayer=1;
102 else if (subType[3]=='P') mLayer=2;
103 }
104 else if (subType[1]=='L')
105 {
106 if (subType[3]=='P') mLayer=1;
107 else if (subType[3]=='C') mLayer=2;
108 }
109
110 if (mLayer != 1 && mLayer !=2) G4cout << " something is wrong - multilayer index is " << mLayer << G4endl;
111
112 int sTgcId = m_muonHelper->BuildsTgcHitId(subType, iPhi, iRing, mLayer,nLayer, iSide);
113 TrackHelper trHelp(aStep->GetTrack());
114
115 m_sTGCSimHitCollection->Emplace(sTgcId,globalTime,position,pdgCode,direction,depositEnergy,
116 trHelp.GenerateParticleLink(m_g4UserEventInfo ? m_g4UserEventInfo->GetEventStore() : nullptr),
117 preStep->GetKineticEnergy(),preposition);
118
119 return true;
120}
#define M_PI
#define ATH_MSG_VERBOSE(x)
#define ATH_MSG_WARNING(x)
double charge(const T &p)
Definition AtlasPID.h:997
AthMessaging(IMessageSvc *msgSvc, const std::string &name)
Constructor.
static AtlasG4EventUserInfo * GetEventUserInfo()
HepMcParticleLink GenerateParticleLink()
Generates a creates new HepMcParticleLink object on the stack based on GetUniqueID(),...
Definition TrackHelper.h:52
unsigned m_baseDepth
basic depth to travel along the G4 history.
sTGCSimHitCollection * m_sTGCSimHitCollection
AtlasG4EventUserInfo * m_g4UserEventInfo
G4bool ProcessHits(G4Step *aStep, G4TouchableHistory *ROhist) override final
sTGCSensitiveDetector(const std::string &name, const std::string &hitCollectionName, unsigned baseDepth)
construction/destruction
const sTgcHitIdHelper * m_muonHelper
void Initialize(G4HCofThisEvent *HCE) override final
member functions
std::string toString(const Translation3D &translation, int precision=4)
GeoPrimitvesToStringConverter.
Amg::Vector3D Hep3VectorToEigen(const CLHEP::Hep3Vector &CLHEPvector)
Converts a CLHEP-based CLHEP::Hep3Vector into an Eigen-based Amg::Vector3D.
Eigen::Matrix< double, 3, 1 > Vector3D
AtlasHitsVector< sTGCSimHit > sTGCSimHitCollection