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

#include <MDTSensitiveDetector.h>

Inheritance diagram for MDTSensitiveDetector:
Collaboration diagram for MDTSensitiveDetector:

Public Member Functions

 MDTSensitiveDetector (const std::string &name, const std::string &hitCollectionName, const unsigned int nTubesMax)
 construction/destruction
void Initialize (G4HCofThisEvent *HCE) override final
 member functions
G4bool ProcessHits (G4Step *aStep, G4TouchableHistory *ROhist) override final

Private Member Functions

 FRIEND_TEST (MDTSensitiveDetectortest, Initialize)
 FRIEND_TEST (MDTSensitiveDetectortest, ProcessHits)
 FRIEND_TEST (MDTSensitiveDetectortest, GetIdentifier)
int GetIdentifier (const G4TouchableHistory *touchHist)

Private Attributes

std::string m_hitCollectionName
 member data
MDTSimHitCollectionm_MDTHitColl {nullptr}
AtlasG4EventUserInfom_g4UserEventInfo {nullptr}
const MdtHitIdHelperm_muonHelper
double m_driftRadius
double m_globalTime
Amg::Vector3D m_localPosition
double m_DEFAULT_TUBE_RADIUS
 radius assigned to radius if radius is invalid

Detailed Description

Author
danie.nosp@m.la.r.nosp@m.ebuzz.nosp@m.i@pv.nosp@m..infn.nosp@m..it
nveld.nosp@m.ik@n.nosp@m.ikhef.nosp@m..nl

Class methods and properties

The method MDTSensitiveDetector::ProcessHits is executed by the G4 kernel each time a charged particle (or a geantino) crosses one of the MDT Sensitive Gas volumes.

Once a G4Step is perfomed by the particle in the sensitive volume (both when the particle leaves the tube or stops in it), Pre and PostStepPositions are tranformed into local coordinates (chamber reference system, with Z along the tube direction and XY tranversal plane) and used to calculate the local direction of the track.

Two cases are given: 1) the particle over-passed the wire: here the drift radius (the impact parameter) is computed as the local direction distance from the wire;
2) the step doesn't allow the particle to over-pass the wire: here the shortest distance to the wire, which is one of the two end-points, is calculated for each step which occurs inside the sensitive volumes and only the closer one is kept; this case includes also the hard scattering case inside the sensitive volume.

Navigating with the touchableHistory method GetHistoryDepth() through the hierarchy of the volumes crossed by the particles, the MDTSensitiveDetector determinates the correct set of geometry parameters to be folded in the Simulation Identifier associated to each hit. The MDT SimIDs are 32-bit unsigned integers, built using the MuonSimEvent/MdtHitIdHelper class which inherits from the MuonHitIdHelper base class.

We describe in the following, how each field of the identifier is retrieved.

1) stationName, stationEta, stationPhi: when a volume is found in the hierarchy, whose name contains the substring "station", the stationName is extracted from the volume's name; stationPhi and stationEta are calculated starting from the volume copy number, assigned by MuonGeoModel.

2) multilayer: when a volume is found in the hierarchy, whose name contains the substring "component", the multilayer is set to 1 or to 2, according to the component number (multilayer=1 if the component is 1, 5 or 8; multilayer=2 if the component is 3, 7, 8, 10, 11 or 14).

3) tubeLayer and tube: when a volume is found in the hierarchy, whose name contains the substring "Drift", tubeLayer and tube are calculated starting from the Drift volume copy number.

notes:

1) this implementation of the MDT Sensitive Detectors can handle different situations: hard scattering of the muon on the sensitive volume (its direction changes), soft secondary particles completely included in the sensitive volume, muon hits masked by secondaries, like delta rays.

2) for each hit, the time of flight (the G4 globalTime), is recorded and associated to the hit.

3) more than none MDT hit can occur in the same tube. The hit selection is done at the level of the digitization procedure.

4) the MDTHit object contains: the SimID, the drift radius and the globalTime.

Definition at line 85 of file MDTSensitiveDetector.h.

Constructor & Destructor Documentation

◆ MDTSensitiveDetector()

MDTSensitiveDetector::MDTSensitiveDetector ( const std::string & name,
const std::string & hitCollectionName,
const unsigned int nTubesMax )

construction/destruction

Definition at line 24 of file MDTSensitiveDetector.cxx.

25 : G4VSensitiveDetector( name )
26 , m_hitCollectionName( hitCollectionName )
27 , m_driftRadius(0.)
28 , m_globalTime(0.)
29 , m_DEFAULT_TUBE_RADIUS( std::numeric_limits<double>::max() )
30{
32}
const MdtHitIdHelper * m_muonHelper
double m_DEFAULT_TUBE_RADIUS
radius assigned to radius if radius is invalid
std::string m_hitCollectionName
member data
static const MdtHitIdHelper * GetHelper(unsigned int nTubes=78)

Member Function Documentation

◆ FRIEND_TEST() [1/3]

MDTSensitiveDetector::FRIEND_TEST ( MDTSensitiveDetectortest ,
GetIdentifier  )
private

◆ FRIEND_TEST() [2/3]

MDTSensitiveDetector::FRIEND_TEST ( MDTSensitiveDetectortest ,
Initialize  )
private

◆ FRIEND_TEST() [3/3]

MDTSensitiveDetector::FRIEND_TEST ( MDTSensitiveDetectortest ,
ProcessHits  )
private

◆ GetIdentifier()

int MDTSensitiveDetector::GetIdentifier ( const G4TouchableHistory * touchHist)
private

Definition at line 145 of file MDTSensitiveDetector.cxx.

146{
147 // attributes of the MDT identifier construction
148 std::string stationName;
149 int stationEta=1;
150 int stationPhi=1;
151 int multilayer=1;
152 int tubeLayer=1;
153 int tube=1;
154
155 bool isAssembly = false;
156 // scan geometry tree to identify the tube
157 for (int i = touchHist->GetHistoryDepth(); i>=0; i--) {
158
159 std::string::size_type npos;
160 std::string::size_type loc1;
161 std::string::size_type loc2;
162 std::string volName = touchHist->GetVolume(i)->GetName();
163
164 // check if this station is an assembly
165 if ((npos = volName.find("av_")) != std::string::npos &&
166 (npos = volName.find("impr_")) != std::string::npos) isAssembly = true;
167
168 // station: name, eta and phi (-> chamber!)
169 if ((npos = volName.find("station")) != std::string::npos && (!isAssembly)) {
170
171 volName.resize(npos-2);
172 int volCopyNo = touchHist->GetVolume(i)->GetCopyNo();
173 volCopyNo=volCopyNo%1000;
174 stationName = volName;
175 stationEta = volCopyNo/100;
176 stationPhi = abs(volCopyNo%100);
177
178 }
179 else if ((npos = volName.find("component")) != std::string::npos && (!isAssembly)) { // multilayer
180
181 int gmID = 0;
182 if ((loc1 = volName.find('[')) != std::string::npos) {
183 if ((loc2 = volName.find(']', loc1+1)) != std::string::npos) {
184 std::istringstream istrvar(volName.substr(loc1+1,loc2-loc1-1));
185 istrvar>>gmID;
186 }
187 }
188 multilayer = gmID;
189
190 } else if ((npos = volName.find("MDT")) != std::string::npos && isAssembly) {
191
192 // vol name for Assembly components are
193 // av_WWW_impr_XXX_Muon::BMSxMDTxx_pv_ZZZ_NAME
194 // where WWW is ass. istance nr.
195 // XXX is comp. imprint nr.
196 // BMSxMDTxx is the name of the comp. log.Vol.
197 // x station sub-type; xx technology subtype
198 // ZZZ is the comp. number of order
199 // NAME is the comp. tag (geoIdentifierTag)
200 // for MDTs NAME is ml[1] or ml[2]
201 // copy numbers for Ass.components are =
202 // CopyNoBase(= geoIdentifierTag of the assembly) + geoIdentifierTag of the component
203 // geoIdentifierTag of the component = Job
204 // geoIdentifierTag of the assembly = (sideC*10000 +
205 // mirsign*1000 + abs(zi)*100 + fi+1)*100000;
206 // mirsign*1000 + abs(zi)*100 + fi+1)*100000;
207 //
208 if ((loc1 = volName.find("Muon::")) != std::string::npos) {
209 stationName = volName.substr(loc1+6,3);
210 }
211
212 int copyNr = touchHist->GetVolume(i)->GetCopyNo();
213 int copyNrBase = int(copyNr/100000);
214 int sideC = int(copyNrBase/10000);
215 int zi = int((copyNrBase%1000)/100);
216 int fi = int(copyNrBase%100);
217 if (sideC == 1) zi = -zi;
218 stationEta = zi;
219 stationPhi = fi;
220
221 int gmID = 0;
222 if ((loc1 = volName.find('[')) != std::string::npos) {
223 if ((loc2 = volName.find(']', loc1+1)) != std::string::npos) {
224 std::istringstream istrvar(volName.substr(loc1+1,loc2-loc1-1));
225 istrvar>>gmID;
226 }
227 }
228 multilayer = gmID;
229 } else if ((npos = volName.find("Drift")) != std::string::npos) { // layer and tube
230 tubeLayer = touchHist->GetVolume(i)->GetCopyNo()/maxNTubesPerLayer;
231 tube = touchHist->GetVolume(i)->GetCopyNo()%maxNTubesPerLayer;
232 }
233 }
234 //construct the hit identifier
235 return m_muonHelper->BuildMdtHitId(stationName, stationPhi, stationEta, multilayer,tubeLayer, tube);
236}
constexpr uint8_t stationPhi
station Phi 1 to 8
unsigned int constexpr sideC
Definition RPDUtils.h:15
@ loc2
generic first and second local coordinate
Definition ParamDefs.h:35
@ loc1
Definition ParamDefs.h:34

◆ Initialize()

void MDTSensitiveDetector::Initialize ( G4HCofThisEvent * HCE)
finaloverride

member functions

Definition at line 35 of file MDTSensitiveDetector.cxx.

36{
37 m_MDTHitColl = nullptr;
38 if (auto* eventInfo = AtlasG4EventUserInfo::GetEventUserInfo()) {
39 m_MDTHitColl = eventInfo->GetHitCollectionMap()->Find<MDTSimHitCollection>(m_hitCollectionName);
40 m_g4UserEventInfo = eventInfo;
41 }
43}
AtlasHitsVector< MDTSimHit > MDTSimHitCollection
static AtlasG4EventUserInfo * GetEventUserInfo()
AtlasG4EventUserInfo * m_g4UserEventInfo
MDTSimHitCollection * m_MDTHitColl

◆ ProcessHits()

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

Definition at line 45 of file MDTSensitiveDetector.cxx.

45 {
46
47 if (!m_MDTHitColl) {
48 G4Exception("MDTSensitiveDetector::ProcessHits", "MDTHitCollectionMissing", FatalException,
49 "Hit collection not initialized; did SetupEvent run?");
50 return false;
51 }
52 G4Track* currentTrack = aStep->GetTrack();
53
54 // MDTs sensitive to charged particle only
55 if (currentTrack->GetDefinition()->GetPDGCharge() == 0.0) {
56 if (currentTrack->GetDefinition()!=G4Geantino::GeantinoDefinition()) return true;
57 else if (currentTrack->GetDefinition()==G4ChargedGeantino::ChargedGeantinoDefinition()) return true;
58 }
59
60 G4VPhysicalVolume* physVolPostStep = aStep->GetPostStepPoint()->GetPhysicalVolume();
61 if (nullptr == physVolPostStep) return true;
62
63 // hit information to be recorded
64 double globalTime;
65 double driftRadius;
66 Amg::Vector3D localPosition;
67 localPosition.setZero();
68
69 // get top transformation
70 const G4AffineTransform trans = currentTrack->GetTouchable()->GetHistory()->GetTopTransform();
71
72 // transform pre and post step positions to local positions
73
74 Amg::Vector3D localVertex1( Amg::Hep3VectorToEigen( trans.TransformPoint(aStep->GetPreStepPoint()->GetPosition()) ) );
75 Amg::Vector3D localVertex2( Amg::Hep3VectorToEigen( trans.TransformPoint(aStep->GetPostStepPoint()->GetPosition()) ) );
76
77 // calculate local direction from begin- and end-point of the step
78 Amg::Vector3D localDirection( (localVertex2 - localVertex1) ); // normalized
79 localDirection.z() = 0.; // look in xy-plane
80
81 // See if particle passed wire by projecting begin- and end-point on the step direction
82 if( (localVertex2.dot(localDirection)) * (localVertex1.dot(localDirection)) < 0 ) { // particle passed wire
83
84 // compute drift radius ( = impact parameter)
85 double Xpos = localVertex1[0];
86 double Ypos = localVertex1[1];
87 double Xdir = localDirection[0];
88 double Ydir = localDirection[1];
89
90 double Alpha = -1*(Xpos*Xdir + Ypos*Ydir)/(Xdir*Xdir + Ydir*Ydir); // localPosition*localDirection
91
92 localPosition = localVertex1 + Alpha*(localVertex2-localVertex1);
93 driftRadius = localPosition.perp();
94 globalTime = aStep->GetPreStepPoint()->GetGlobalTime(); // take pre step time
95
96 }else{ // particle didn't pass wire, one of the end-points is the shortest distance to the wire
97
98 // calculate which of the two end-points is closer to the wire
99 double dist1 = localVertex1.perp();
100 double dist2 = localVertex2.perp();
101 if( dist1 < dist2 ){ // first is closer
102 driftRadius = dist1;
103 localPosition = localVertex1;
104 globalTime = aStep->GetPreStepPoint()->GetGlobalTime();
105 }else{ // second is closer
106 driftRadius = dist2;
107 localPosition = localVertex2;
108 globalTime = aStep->GetPostStepPoint()->GetGlobalTime();
109 }
110 }
111
112 if( driftRadius < m_driftRadius ){ // check if current step came closer to the wire than the previous ones
114 m_globalTime = globalTime;
115 m_localPosition = localPosition;
116 }
117
118 // check if particle left tube or stopped in tube
119 G4String namePreStepMat = aStep->GetPreStepPoint()->GetMaterial()->GetName();
120 G4String namePostStepMat = aStep->GetPostStepPoint()->GetMaterial()->GetName();
121 G4String nameSD = aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetSensitiveDetector()->GetName();
122 // G4int trackid = aStep->GetTrack()->GetTrackID();
123 // see if we were in the sensitive volume and left it, or the particle was killed
124 if( ((nameSD) && (namePreStepMat != namePostStepMat)) || (currentTrack->GetTrackStatus() == fStopAndKill)){
125
126 // get identifier
127 const G4TouchableHistory* touchHist = static_cast<const G4TouchableHistory*>(aStep->GetPreStepPoint()->GetTouchable());
128 int MDTid = GetIdentifier(touchHist);
129
130 TrackHelper trHelp(aStep->GetTrack());
131
132 // construct new mdt hit
134 trHelp.GenerateParticleLink(m_g4UserEventInfo ? m_g4UserEventInfo->GetEventStore() : nullptr),
135 aStep->GetStepLength(),
136 aStep->GetTotalEnergyDeposit(),
137 currentTrack->GetDefinition()->GetPDGEncoding(),
138 aStep->GetPreStepPoint()->GetKineticEnergy());
139
140 m_driftRadius = m_DEFAULT_TUBE_RADIUS; // reset start value of driftRadius
141 }
142 return true;
143}
int GetIdentifier(const G4TouchableHistory *touchHist)
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
@ driftRadius
trt, straws
Definition ParamDefs.h:53

Member Data Documentation

◆ m_DEFAULT_TUBE_RADIUS

double MDTSensitiveDetector::m_DEFAULT_TUBE_RADIUS
private

radius assigned to radius if radius is invalid

Definition at line 112 of file MDTSensitiveDetector.h.

◆ m_driftRadius

double MDTSensitiveDetector::m_driftRadius
private

Definition at line 107 of file MDTSensitiveDetector.h.

◆ m_g4UserEventInfo

AtlasG4EventUserInfo* MDTSensitiveDetector::m_g4UserEventInfo {nullptr}
private

Definition at line 104 of file MDTSensitiveDetector.h.

104{nullptr};

◆ m_globalTime

double MDTSensitiveDetector::m_globalTime
private

Definition at line 108 of file MDTSensitiveDetector.h.

◆ m_hitCollectionName

std::string MDTSensitiveDetector::m_hitCollectionName
private

member data

Definition at line 102 of file MDTSensitiveDetector.h.

◆ m_localPosition

Amg::Vector3D MDTSensitiveDetector::m_localPosition
private

Definition at line 109 of file MDTSensitiveDetector.h.

◆ m_MDTHitColl

MDTSimHitCollection* MDTSensitiveDetector::m_MDTHitColl {nullptr}
private

Definition at line 103 of file MDTSensitiveDetector.h.

103{nullptr};

◆ m_muonHelper

const MdtHitIdHelper* MDTSensitiveDetector::m_muonHelper
private

Definition at line 105 of file MDTSensitiveDetector.h.


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