ATLAS Offline Software
Loading...
Searching...
No Matches
GenericMuonSensitiveDetector.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2019 CERN for the benefit of the ATLAS collaboration
3*/
4
8
9
10#include "G4Exception.hh"
11#include "G4Track.hh"
12
13#include <string>
14
16
17// construction/destruction
18GenericMuonSensitiveDetector::GenericMuonSensitiveDetector(const std::string& name, const std::string& hitCollectionName)
19 : G4VSensitiveDetector( name )
20 , m_hitCollectionName( hitCollectionName )
21{
22 G4cout << " creating a GenericMuonSensitiveDetector: "<<name << G4endl;
23}
24
25// Implemenation of memebr functions
27{
29 if (auto* eventInfo = AtlasG4EventUserInfo::GetEventUserInfo()) {
30 m_GenericMuonHitCollection = eventInfo->GetHitCollectionMap()->Find<GenericMuonSimHitCollection>(m_hitCollectionName);
31 m_g4UserEventInfo = eventInfo;
32 }
33}
34
35G4bool GenericMuonSensitiveDetector::ProcessHits(G4Step* aStep,G4TouchableHistory* /*ROHist*/)
36{
38 G4Exception("GenericMuonSensitiveDetector::ProcessHits", "GenericMuonHitCollectionMissing", FatalException,
39 "Hit collection not initialized; did SetupEvent run?");
40 return false;
41 }
42 G4cout << "Hit in a sensitive layer!!!!! " << G4endl;
43 G4Track* currentTrack = aStep->GetTrack();
44 const G4AffineTransform trans = currentTrack->GetTouchable()->GetHistory()->GetTopTransform(); // from global to local
45 G4StepPoint* postStep=aStep->GetPostStepPoint();
46 G4StepPoint* preStep=aStep->GetPreStepPoint();
47 const G4Step* post_Step=aStep->GetTrack()->GetStep();
48
49 Amg::Vector3D position = Amg::Hep3VectorToEigen( postStep->GetPosition() );
50 Amg::Vector3D preposition = Amg::Hep3VectorToEigen( preStep->GetPosition() );
51
52 Amg::Vector3D local_position = Amg::Hep3VectorToEigen( trans.TransformPoint( postStep->GetPosition() ) ) ;
53 Amg::Vector3D local_preposition = Amg::Hep3VectorToEigen( trans.TransformPoint( preStep->GetPosition() ) );
54
55 int pdgCode=currentTrack->GetDefinition()->GetPDGEncoding();
56
57 float globalTime=postStep->GetGlobalTime();
58 float globalpreTime=preStep->GetGlobalTime();
59 float eKin=postStep->GetKineticEnergy();
60
61 Amg::Vector3D direction = Amg::Hep3VectorToEigen( postStep->GetMomentumDirection() );
62 float depositEnergy=post_Step->GetTotalEnergyDeposit();
63 float StepLength=post_Step->GetStepLength();
64
65 TrackHelper trHelp(aStep->GetTrack());
66
67 //G4cout << aHit->print() << G4endl;
68 m_GenericMuonHitCollection->Emplace( 0 /* HitID id generic*/,globalTime,globalpreTime,position,local_position,preposition,local_preposition,pdgCode,eKin,direction,depositEnergy,StepLength,
69 trHelp.GenerateParticleLink(m_g4UserEventInfo ? m_g4UserEventInfo->GetEventStore() : nullptr));
70
71 return true;
72}
AtlasHitsVector< GenericMuonSimHit > GenericMuonSimHitCollection
static AtlasG4EventUserInfo * GetEventUserInfo()
G4bool ProcessHits(G4Step *aStep, G4TouchableHistory *ROhist) override final
GenericMuonSensitiveDetector(const std::string &name, const std::string &hitCollectionName)
construction/destruction
void Initialize(G4HCofThisEvent *HCE) override final
member functions
GenericMuonSimHitCollection * m_GenericMuonHitCollection
HepMcParticleLink GenerateParticleLink()
Generates a creates new HepMcParticleLink object on the stack based on GetUniqueID(),...
Definition TrackHelper.h:52
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