ATLAS Offline Software
Public Member Functions | Private Member Functions | Private Attributes | List of all members
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)
 construction/destruction More...
 
 ~sTGCSensitiveDetector ()
 
void Initialize (G4HCofThisEvent *HCE) override final
 member functions More...
 
G4bool ProcessHits (G4Step *aStep, G4TouchableHistory *ROhist) override final
 

Private Member Functions

 FRIEND_TEST (sTGCSensitiveDetectortest, Initialize)
 
 FRIEND_TEST (sTGCSensitiveDetectortest, ProcessHits)
 

Private Attributes

SG::WriteHandle< sTGCSimHitCollectionm_sTGCSimHitCollection
 
const sTgcHitIdHelperm_muonHelper
 

Detailed Description

Definition at line 15 of file sTGCSensitiveDetector.h.

Constructor & Destructor Documentation

◆ sTGCSensitiveDetector()

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

construction/destruction

Definition at line 18 of file sTGCSensitiveDetector.cxx.

19  : G4VSensitiveDetector( name )
20  , m_sTGCSimHitCollection( hitCollectionName )
21 {
23  //m_muonHelper->PrintFields();
24 }

◆ ~sTGCSensitiveDetector()

sTGCSensitiveDetector::~sTGCSensitiveDetector ( )
inline

Definition at line 22 of file sTGCSensitiveDetector.h.

22 {}

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 27 of file sTGCSensitiveDetector.cxx.

28 {
29  if (!m_sTGCSimHitCollection.isValid()) m_sTGCSimHitCollection = std::make_unique<sTGCSimHitCollection>();
30 }

◆ ProcessHits()

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

Definition at line 32 of file sTGCSensitiveDetector.cxx.

33 {
34  G4Track* currentTrack = aStep->GetTrack();
35  int charge=currentTrack->GetDefinition()->GetPDGCharge();
36 
37  bool geantinoHit = (currentTrack->GetDefinition()==G4Geantino::GeantinoDefinition()) ||
38  (currentTrack->GetDefinition()==G4ChargedGeantino::ChargedGeantinoDefinition());
39 
40  if (!charge && (!geantinoHit)) return false;
41  // G4cout << "\t\t sTGCSD: Hit in a sensitive layer!!!!! " << G4endl;
42  G4StepPoint* postStep=aStep->GetPostStepPoint();
43  const G4Step* post_Step=aStep->GetTrack()->GetStep();
44 
45  Amg::Vector3D position = Amg::Hep3VectorToEigen(postStep->GetPosition());
46 
47  G4StepPoint* preStep = aStep->GetPreStepPoint();
48  Amg::Vector3D preposition = Amg::Hep3VectorToEigen(preStep->GetPosition());
49 
50  int pdgCode=currentTrack->GetDefinition()->GetPDGEncoding();
51 
52  float globalTime=postStep->GetGlobalTime();
53 
54  Amg::Vector3D direction= Amg::Hep3VectorToEigen( postStep->GetMomentumDirection() );
55  float depositEnergy=post_Step->GetTotalEnergyDeposit();
56 
57  if (depositEnergy<0.0001 && (!geantinoHit)) return false;
58 
59  const G4TouchableHistory* touchHist = static_cast<const G4TouchableHistory*>(aStep->GetPreStepPoint()->GetTouchable());
60 
61  // int iDepth=touchHist->GetHistoryDepth();
62  // G4cout << "\t\t\t\t Touchable history dump "<<G4endl;
63  int nLayer=touchHist->GetVolume(0)->GetCopyNo();
64  std::string chName=touchHist->GetVolume(1)->GetLogicalVolume()->GetName();
65  //G4cout << "sTGCSensitiveDetector name: "<<chName<<G4endl;
66  std::string subType=chName.substr(chName.find('-')+1);
67  //G4cout << "\t\t sType: "<<subType);
68  if (subType[0]!='T'&&subType[0]!='Q' ) G4cout << " something is wrong, this is no sTGC!"<<G4endl;
69  std::string temp(&subType[2]);
70  std::istringstream is(temp);
71  int iRing;
72  is>>iRing;
73  // identifiers have eta naming 0-2, eta encoded in subtype is 1-3
74  iRing--;
75  // double phiDiff=2*M_PI;
76 
77  G4ThreeVector posH=postStep->GetPosition(); //posH is equivalent to position - eigen not used to avoid additional dependence on EventPrimitives
78  if (subType[1]=='L') posH.rotateZ(M_PI/8.);
79  double phiHit=posH.phi();
80  if(phiHit<=0) phiHit+=2.*M_PI;
81  int iPhi=1+int(phiHit/(M_PI/4.));
82  iPhi*=2;
83  if (subType[1]=='L') iPhi-=1;
84 
85  int iSide=1;
86  if (position.z()<0) iSide=-1;
87 
88  int mLayer=0;
89  if (subType[1]=='S')
90  {
91  if (subType[3]=='C') mLayer=1;
92  else if (subType[3]=='P') mLayer=2;
93  }
94  else if (subType[1]=='L')
95  {
96  if (subType[3]=='P') mLayer=1;
97  else if (subType[3]=='C') mLayer=2;
98  }
99 
100  if (mLayer != 1 && mLayer !=2) G4cout << " something is wrong - multilayer index is " << mLayer << G4endl;
101 
102  int sTgcId = m_muonHelper->BuildsTgcHitId(subType, iPhi, iRing, mLayer,nLayer, iSide);
103  TrackHelper trHelp(aStep->GetTrack());
104  m_sTGCSimHitCollection->Emplace(sTgcId,globalTime,position,pdgCode,direction,depositEnergy,
105  trHelp.GenerateParticleLink(),
106  preStep->GetKineticEnergy(),preposition);
107 
108  return true;
109 }

Member Data Documentation

◆ m_muonHelper

const sTgcHitIdHelper* sTGCSensitiveDetector::m_muonHelper
private

Definition at line 31 of file sTGCSensitiveDetector.h.

◆ m_sTGCSimHitCollection

SG::WriteHandle<sTGCSimHitCollection> sTGCSensitiveDetector::m_sTGCSimHitCollection
private

Definition at line 30 of file sTGCSensitiveDetector.h.


The documentation for this class was generated from the following files:
CaloCellPos2Ntuple.int
int
Definition: CaloCellPos2Ntuple.py:24
sTgcHitIdHelper::GetHelper
static const sTgcHitIdHelper * GetHelper()
Definition: sTgcHitIdHelper.cxx:24
M_PI
#define M_PI
Definition: ActiveFraction.h:11
sTGCSensitiveDetector::m_muonHelper
const sTgcHitIdHelper * m_muonHelper
Definition: sTGCSensitiveDetector.h:31
Amg::Hep3VectorToEigen
Amg::Vector3D Hep3VectorToEigen(const CLHEP::Hep3Vector &CLHEPvector)
Converts a CLHEP-based CLHEP::Hep3Vector into an Eigen-based Amg::Vector3D.
Definition: CLHEPtoEigenConverter.h:137
sTgcHitIdHelper::BuildsTgcHitId
int BuildsTgcHitId(const std::string &, const int, const int, const int, const int, const int) const
Definition: sTgcHitIdHelper.cxx:89
TrackHelper
Definition: TrackHelper.h:14
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:195
charge
double charge(const T &p)
Definition: AtlasPID.h:494
Amg::Vector3D
Eigen::Matrix< double, 3, 1 > Vector3D
Definition: GeoPrimitives.h:47
Trk::iPhi
@ iPhi
Definition: ParamDefs.h:53
sTGCSensitiveDetector::m_sTGCSimHitCollection
SG::WriteHandle< sTGCSimHitCollection > m_sTGCSimHitCollection
Definition: sTGCSensitiveDetector.h:30