ATLAS Offline Software
Loading...
Searching...
No Matches
sTGCSensitiveDetector Class Reference

#include <sTGCSensitiveDetector.h>

Inheritance diagram for sTGCSensitiveDetector:
Collaboration diagram for sTGCSensitiveDetector:

Public Member Functions

 sTGCSensitiveDetector (const std::string &name, const std::string &hitCollectionName, unsigned baseDepth)
 construction/destruction
void Initialize (G4HCofThisEvent *HCE) override final
 member functions
G4bool ProcessHits (G4Step *aStep, G4TouchableHistory *ROhist) override final
bool msgLvl (const MSG::Level lvl) const
 Test the output level.
MsgStream & msg () const
 The standard message stream.
MsgStream & msg (const MSG::Level lvl) const
 The standard message stream.
void setLevel (MSG::Level lvl)
 Change the current logging level.

Private Member Functions

 FRIEND_TEST (sTGCSensitiveDetectortest, Initialize)
 FRIEND_TEST (sTGCSensitiveDetectortest, ProcessHits)
void initMessaging () const
 Initialize our message level and MessageSvc.

Private Attributes

std::string m_hitCollectionName
sTGCSimHitCollectionm_sTGCSimHitCollection {nullptr}
AtlasG4EventUserInfom_g4UserEventInfo {nullptr}
const sTgcHitIdHelperm_muonHelper {sTgcHitIdHelper::GetHelper()}
unsigned m_baseDepth {0}
 basic depth to travel along the G4 history.
std::string m_nm
 Message source name.
boost::thread_specific_ptr< MsgStream > m_msg_tls
 MsgStream instance (a std::cout like with print-out levels)
std::atomic< IMessageSvc * > m_imsg { nullptr }
 MessageSvc pointer.
std::atomic< MSG::Level > m_lvl { MSG::NIL }
 Current logging level.
std::atomic_flag m_initialized ATLAS_THREAD_SAFE = ATOMIC_FLAG_INIT
 Messaging initialized (initMessaging)

Detailed Description

Definition at line 18 of file sTGCSensitiveDetector.h.

Constructor & Destructor Documentation

◆ sTGCSensitiveDetector()

sTGCSensitiveDetector::sTGCSensitiveDetector ( const std::string & name,
const std::string & hitCollectionName,
unsigned baseDepth )

construction/destruction

Definition at line 21 of file sTGCSensitiveDetector.cxx.

23 :
24 G4VSensitiveDetector( name ),
26 m_hitCollectionName( hitCollectionName ),
27 m_baseDepth{baseDepth} {}
AthMessaging(IMessageSvc *msgSvc, const std::string &name)
Constructor.
unsigned m_baseDepth
basic depth to travel along the G4 history.

Member Function Documentation

◆ FRIEND_TEST() [1/2]

sTGCSensitiveDetector::FRIEND_TEST ( sTGCSensitiveDetectortest ,
Initialize  )
private

◆ FRIEND_TEST() [2/2]

sTGCSensitiveDetector::FRIEND_TEST ( sTGCSensitiveDetectortest ,
ProcessHits  )
private

◆ Initialize()

void sTGCSensitiveDetector::Initialize ( G4HCofThisEvent * HCE)
finaloverride

member functions

Definition at line 30 of file sTGCSensitiveDetector.cxx.

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}
static AtlasG4EventUserInfo * GetEventUserInfo()
sTGCSimHitCollection * m_sTGCSimHitCollection
AtlasG4EventUserInfo * m_g4UserEventInfo
AtlasHitsVector< sTGCSimHit > sTGCSimHitCollection

◆ initMessaging()

void AthMessaging::initMessaging ( ) const
privateinherited

Initialize our message level and MessageSvc.

This method should only be called once.

Definition at line 39 of file AthMessaging.cxx.

40{
42 // If user did not set an explicit level, set a default
43 if (m_lvl == MSG::NIL) {
44 m_lvl = m_imsg ?
45 static_cast<MSG::Level>( m_imsg.load()->outputLevel(m_nm) ) :
46 MSG::INFO;
47 }
48}
std::string m_nm
Message source name.
std::atomic< IMessageSvc * > m_imsg
MessageSvc pointer.
std::atomic< MSG::Level > m_lvl
Current logging level.
IMessageSvc * getMessageSvc(bool quiet=false)

◆ msg() [1/2]

MsgStream & AthMessaging::msg ( ) const
inlineinherited

The standard message stream.

Returns a reference to the default message stream May not be invoked before sysInitialize() has been invoked.

Definition at line 163 of file AthMessaging.h.

164{
165 MsgStream* ms = m_msg_tls.get();
166 if (!ms) {
167 if (!m_initialized.test_and_set()) initMessaging();
168 ms = new MsgStream(m_imsg,m_nm);
169 m_msg_tls.reset( ms );
170 }
171
172 ms->setLevel (m_lvl);
173 return *ms;
174}
boost::thread_specific_ptr< MsgStream > m_msg_tls
MsgStream instance (a std::cout like with print-out levels)
void initMessaging() const
Initialize our message level and MessageSvc.

◆ msg() [2/2]

MsgStream & AthMessaging::msg ( const MSG::Level lvl) const
inlineinherited

The standard message stream.

Returns a reference to the default message stream May not be invoked before sysInitialize() has been invoked.

Definition at line 178 of file AthMessaging.h.

179{ return msg() << lvl; }
MsgStream & msg() const
The standard message stream.

◆ msgLvl()

bool AthMessaging::msgLvl ( const MSG::Level lvl) const
inlineinherited

Test the output level.

Parameters
lvlThe message level to test against
Returns
boolean Indicating if messages at given level will be printed
Return values
trueMessages at level "lvl" will be printed

Definition at line 151 of file AthMessaging.h.

152{
153 if (m_lvl <= lvl) {
154 msg() << lvl;
155 return true;
156 } else {
157 return false;
158 }
159}

◆ ProcessHits()

G4bool sTGCSensitiveDetector::ProcessHits ( G4Step * aStep,
G4TouchableHistory * ROhist )
finaloverride

Definition at line 39 of file sTGCSensitiveDetector.cxx.

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
const sTgcHitIdHelper * m_muonHelper
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
const std::string & chName(ChIndex index)
convert ChIndex into a string
@ iPhi
Definition ParamDefs.h:47

◆ setLevel()

void AthMessaging::setLevel ( MSG::Level lvl)
inherited

Change the current logging level.

Use this rather than msg().setLevel() for proper operation with MT.

Definition at line 28 of file AthMessaging.cxx.

29{
30 m_lvl = lvl;
31}

Member Data Documentation

◆ ATLAS_THREAD_SAFE

std::atomic_flag m_initialized AthMessaging::ATLAS_THREAD_SAFE = ATOMIC_FLAG_INIT
mutableprivateinherited

Messaging initialized (initMessaging)

Definition at line 141 of file AthMessaging.h.

◆ m_baseDepth

unsigned sTGCSensitiveDetector::m_baseDepth {0}
private

basic depth to travel along the G4 history.

For jobs run with the legacy geometry database, it's zero. Otherwise, in the new sqlite workflow it's 1

Definition at line 41 of file sTGCSensitiveDetector.h.

41{0};

◆ m_g4UserEventInfo

AtlasG4EventUserInfo* sTGCSensitiveDetector::m_g4UserEventInfo {nullptr}
private

Definition at line 37 of file sTGCSensitiveDetector.h.

37{nullptr};

◆ m_hitCollectionName

std::string sTGCSensitiveDetector::m_hitCollectionName
private

Definition at line 35 of file sTGCSensitiveDetector.h.

◆ m_imsg

std::atomic<IMessageSvc*> AthMessaging::m_imsg { nullptr }
mutableprivateinherited

MessageSvc pointer.

Definition at line 135 of file AthMessaging.h.

135{ nullptr };

◆ m_lvl

std::atomic<MSG::Level> AthMessaging::m_lvl { MSG::NIL }
mutableprivateinherited

Current logging level.

Definition at line 138 of file AthMessaging.h.

138{ MSG::NIL };

◆ m_msg_tls

boost::thread_specific_ptr<MsgStream> AthMessaging::m_msg_tls
mutableprivateinherited

MsgStream instance (a std::cout like with print-out levels)

Definition at line 132 of file AthMessaging.h.

◆ m_muonHelper

const sTgcHitIdHelper* sTGCSensitiveDetector::m_muonHelper {sTgcHitIdHelper::GetHelper()}
private

Definition at line 38 of file sTGCSensitiveDetector.h.

static const sTgcHitIdHelper * GetHelper()

◆ m_nm

std::string AthMessaging::m_nm
privateinherited

Message source name.

Definition at line 129 of file AthMessaging.h.

◆ m_sTGCSimHitCollection

sTGCSimHitCollection* sTGCSensitiveDetector::m_sTGCSimHitCollection {nullptr}
private

Definition at line 36 of file sTGCSensitiveDetector.h.

36{nullptr};

The documentation for this class was generated from the following files: