ATLAS Offline Software
Loading...
Searching...
No Matches
MicromegasSensitiveDetector.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2022 CERN for the benefit of the ATLAS collaboration
3*/
4
9#include "G4Exception.hh"
10#include "G4Geantino.hh"
11#include "G4ChargedGeantino.hh"
12
13#include "G4Track.hh"
14
16
17#include <string>
18
19// construction/destruction
20MicromegasSensitiveDetector::MicromegasSensitiveDetector(const std::string& name, const std::string& hitCollectionName)
21 : G4VSensitiveDetector( name )
22 , m_hitCollectionName( hitCollectionName )
23{
25 //m_muonHelper->PrintFields();
26}
27
28// Implemenation of memebr functions
30{
31 m_MMSimHitCollection = nullptr;
32 if (auto* eventInfo = AtlasG4EventUserInfo::GetEventUserInfo()) {
33 m_MMSimHitCollection = eventInfo->GetHitCollectionMap()->Find<MMSimHitCollection>(m_hitCollectionName);
34 m_g4UserEventInfo = eventInfo;
35 }
36}
37
38G4bool MicromegasSensitiveDetector::ProcessHits(G4Step* aStep,G4TouchableHistory* /*ROHist*/)
39{
41 G4Exception("MicromegasSensitiveDetector::ProcessHits", "MicromegasHitCollectionMissing", FatalException,
42 "Hit collection not initialized; did SetupEvent run?");
43 return false;
44 }
45 G4Track* currentTrack = aStep->GetTrack();
46 int charge=currentTrack->GetDefinition()->GetPDGCharge();
47 bool geantinoHit = (currentTrack->GetDefinition()==G4Geantino::GeantinoDefinition()) ||
48 (currentTrack->GetDefinition()==G4ChargedGeantino::ChargedGeantinoDefinition());
49
50 if (!charge && (!geantinoHit)) return false;
51 // G4cout << "\t\t MicromegasSD: Hit in a sensitive layer!!!!! " << G4endl;
52 G4StepPoint* postStep=aStep->GetPostStepPoint();
53 const G4Step* post_Step=aStep->GetTrack()->GetStep();
54
55 Amg::Vector3D position = Amg::Hep3VectorToEigen( postStep->GetPosition() );
56
57 int pdgCode=currentTrack->GetDefinition()->GetPDGEncoding();
58
59 float globalTime=postStep->GetGlobalTime();
60 float eKin=postStep->GetKineticEnergy();
61 if (eKin<= 0. && (!geantinoHit)) return false;
62
63 Amg::Vector3D direction = Amg::Hep3VectorToEigen( postStep->GetMomentumDirection() );
64 float depositEnergy=post_Step->GetTotalEnergyDeposit();
65
66 if (depositEnergy<0.0001 && (!geantinoHit)) return false;
67
68 const G4TouchableHistory* touchHist = static_cast<const G4TouchableHistory*>(aStep->GetPreStepPoint()->GetTouchable());
69
70 // int iDepth=touchHist->GetHistoryDepth();
71 // G4cout << "\t\t\t\t Touchable history dump " << G4endl;
72 int nLayer=touchHist->GetVolume(0)->GetCopyNo();
73 std::string chName=touchHist->GetVolume(1)->GetLogicalVolume()->GetName();
74 std::string subType=chName.substr(chName.find('-')+1);
75 if (subType[0]!='M') G4cout << " something is wrong, this is no Micromegas!" << G4endl;
76 std::string temp(&subType[1]);
77 std::istringstream is(temp);
78 int iRing;
79 is>>iRing;
80 // identifiers have eta naming 0-1, eta encoded in subtype is 1-2
81 iRing--;
82 // double phiDiff=2*M_PI;
83
84 G4ThreeVector posH=postStep->GetPosition(); //posH is equivalent to position - eigen not used to avoid additional dependence on EventPrimitives
85 if (subType[2]=='L') posH.rotateZ(M_PI/8.);
86 double phiHit=posH.phi();
87 if(phiHit<=0) phiHit+=2.*M_PI;
88 int iPhi=1+int(phiHit/(M_PI/4.));
89 iPhi*=2;
90 if (subType[2]=='L') iPhi-=1;
91
92 int iSide=1;
93 if (position.z()<0) iSide=-1;
94
95 int mLayer= atoi((subType.substr(3,1)).c_str());
96 if (mLayer != 1 && mLayer !=2) G4cout << " something is wrong - multilayer index is " << mLayer << G4endl;
97
98 // G4cout << "\t\t Chamber "<<chName<<" subType "<<subType<<" layer nr. "<<nLayer<<" ring "<<iRing<<" sector "<<iPhi<<" side "<<iSide << G4endl;
99 int MmId = m_muonHelper->BuildMicromegasHitId(subType, iPhi, iRing, mLayer,nLayer, iSide);
100
101 TrackHelper trHelp(aStep->GetTrack());
102
103 m_MMSimHitCollection->Emplace(MmId, globalTime,position,pdgCode,eKin,direction,depositEnergy,
104 trHelp.GenerateParticleLink(m_g4UserEventInfo ? m_g4UserEventInfo->GetEventStore() : nullptr));
105
106 // G4cout << "MMs "<<m_muonHelper->GetStationName(MmId)
107 // << " "<<m_muonHelper->GetFieldValue("PhiSector")
108 // << " "<<m_muonHelper->GetFieldValue("ZSector")
109 // << " "<<m_muonHelper->GetFieldValue("MultiLayer")
110 // << " "<<m_muonHelper->GetFieldValue("Layer")
111 // << " "<<m_muonHelper->GetFieldValue("Side") << G4endl;
112
113 //G4cout << m_muonHelper->GetStationName(MmId)<<" "<<aHit->print() << G4endl;
114 // G4cout << aHit->print() << G4endl;
115
116 return true;
117}
#define M_PI
double charge(const T &p)
Definition AtlasPID.h:997
AtlasHitsVector< MMSimHit > MMSimHitCollection
static AtlasG4EventUserInfo * GetEventUserInfo()
static const MicromegasHitIdHelper * GetHelper()
void Initialize(G4HCofThisEvent *HCE) override final
member functions
const MicromegasHitIdHelper * m_muonHelper
MicromegasSensitiveDetector(const std::string &name, const std::string &hitCollectionName)
construction/destruction
G4bool ProcessHits(G4Step *aStep, G4TouchableHistory *ROhist) override final
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