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

#include <MDTSensitiveDetectorCosmics.h>

Inheritance diagram for MDTSensitiveDetectorCosmics:
Collaboration diagram for MDTSensitiveDetectorCosmics:

Public Member Functions

 MDTSensitiveDetectorCosmics (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 (MDTSensitiveDetectorCosmicstest, Initialize)
 FRIEND_TEST (MDTSensitiveDetectorCosmicstest, ProcessHits)
 FRIEND_TEST (MDTSensitiveDetectorCosmicstest, GetIdentifier)
int GetIdentifier (const G4TouchableHistory *touchHist)

Private Attributes

Amg::Vector3D m_mom
double m_momMag
Amg::Vector3D m_vertex
Amg::Vector3D m_currVertex
Amg::Vector3D m_globH
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 MDTSensitiveDetectorCosmics::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 MDTSensitiveDetectorCosmics 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.

The modified version (Daniela Rebuzzi, 16/12/2005) calculates the tof for the cosmics, using as t0 the time to the closest approach of the track to (0,0,0). This implementation works also for slow particles, since the discrimination physics event/cosmics is done relying on the G4Track vertex (not on the tof).

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 90 of file MDTSensitiveDetectorCosmics.h.

Constructor & Destructor Documentation

◆ MDTSensitiveDetectorCosmics()

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

construction/destruction

Definition at line 29 of file MDTSensitiveDetectorCosmics.cxx.

30 : G4VSensitiveDetector( name )
31 , m_momMag(0.)
32 , m_hitCollectionName( hitCollectionName )
33 , m_driftRadius(0.)
34 , m_globalTime(0.)
35 , m_DEFAULT_TUBE_RADIUS( std::numeric_limits<double>::max() )
36{
38}
double m_DEFAULT_TUBE_RADIUS
radius assigned to radius if radius is invalid
static const MdtHitIdHelper * GetHelper(unsigned int nTubes=78)

Member Function Documentation

◆ FRIEND_TEST() [1/3]

MDTSensitiveDetectorCosmics::FRIEND_TEST ( MDTSensitiveDetectorCosmicstest ,
GetIdentifier  )
private

◆ FRIEND_TEST() [2/3]

MDTSensitiveDetectorCosmics::FRIEND_TEST ( MDTSensitiveDetectorCosmicstest ,
Initialize  )
private

◆ FRIEND_TEST() [3/3]

MDTSensitiveDetectorCosmics::FRIEND_TEST ( MDTSensitiveDetectorCosmicstest ,
ProcessHits  )
private

◆ GetIdentifier()

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

Definition at line 190 of file MDTSensitiveDetectorCosmics.cxx.

191{
192 // attributes of the MDT identifier construction
193 std::string stationName;
194 int stationEta=0;
195 int stationPhi=0;
196 int multilayer=0;
197 int tubeLayer=0;
198 int tube=0;
199
200 bool isAssembly = false;
201 // scan geometry tree to identify the tube
202 for (int i = touchHist->GetHistoryDepth(); i>=0; i--) {
203
204 std::string::size_type npos;
205 std::string::size_type loc1;
206 std::string::size_type loc2;
207 std::string volName = touchHist->GetVolume(i)->GetName();
208
209 // check if this station is an assembly
210 if ((npos = volName.find("av_")) != std::string::npos &&
211 (npos = volName.find("impr_")) != std::string::npos) isAssembly = true;
212
213 // station: name, eta and phi (-> chamber!)
214 if ((npos = volName.find("station")) != std::string::npos && (!isAssembly)) {
215
216 volName.resize(npos-2);
217 int volCopyNo = touchHist->GetVolume(i)->GetCopyNo();
218 volCopyNo=volCopyNo%1000;
219 stationName = volName;
220 stationEta = volCopyNo/100;
221 stationPhi = abs(volCopyNo%100);
222
223 }
224 else if ((npos = volName.find("component")) != std::string::npos && (!isAssembly)) { // multilayer
225
226 int gmID = 0;
227 if ((loc1 = volName.find('[')) != std::string::npos) {
228 if ((loc2 = volName.find(']', loc1+1)) != std::string::npos) {
229 std::istringstream istrvar(volName.substr(loc1+1,loc2-loc1-1));
230 istrvar>>gmID;
231 }
232 }
233 multilayer = gmID;
234
235 } else if ((npos = volName.find("MDT")) != std::string::npos && isAssembly) {
236
237 // vol name for Assembly components are
238 // av_WWW_impr_XXX_Muon::BMSxMDTxx_pv_ZZZ_NAME
239 // where WWW is ass. istance nr.
240 // XXX is comp. imprint nr.
241 // BMSxMDTxx is the name of the comp. log.Vol.
242 // x station sub-type; xx technology subtype
243 // ZZZ is the comp. number of order
244 // NAME is the comp. tag (geoIdentifierTag)
245 // for MDTs NAME is ml[1] or ml[2]
246 // copy numbers for Ass.components are =
247 // CopyNoBase(= geoIdentifierTag of the assembly) + geoIdentifierTag of the component
248 // geoIdentifierTag of the component = Job
249 // geoIdentifierTag of the assembly = (sideC*10000 +
250 // mirsign*1000 + abs(zi)*100 + fi+1)*100000;
251 // mirsign*1000 + abs(zi)*100 + fi+1)*100000;
252 //
253 if ((loc1 = volName.find("Muon::")) != std::string::npos) {
254 stationName = volName.substr(loc1+6,3);
255 }
256
257 int copyNr = touchHist->GetVolume(i)->GetCopyNo();
258 int copyNrBase = int(copyNr/100000);
259 int sideC = int(copyNrBase/10000);
260 int zi = int((copyNrBase%1000)/100);
261 int fi = int(copyNrBase%100);
262 if (sideC == 1) zi = -zi;
263 stationEta = zi;
264 stationPhi = fi;
265
266 int gmID = 0;
267 if ((loc1 = volName.find('[')) != std::string::npos) {
268 if ((loc2 = volName.find(']', loc1+1)) != std::string::npos) {
269 std::istringstream istrvar(volName.substr(loc1+1,loc2-loc1-1));
270 istrvar>>gmID;
271 }
272 }
273 multilayer = gmID;
274 } else if ((npos = volName.find("Drift")) != std::string::npos) { // layer and tube
275 tubeLayer = touchHist->GetVolume(i)->GetCopyNo()/maxNTubesPerLayer;
276 tube = touchHist->GetVolume(i)->GetCopyNo()%maxNTubesPerLayer;
277 }
278 }
279 //construct the hit identifier
280 return m_muonHelper->BuildMdtHitId(stationName, stationPhi, stationEta, multilayer,tubeLayer, tube);
281}
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 MDTSensitiveDetectorCosmics::Initialize ( G4HCofThisEvent * HCE)
finaloverride

member functions

Definition at line 41 of file MDTSensitiveDetectorCosmics.cxx.

42{
43 m_MDTHitColl = nullptr;
44 if (auto* eventInfo = AtlasG4EventUserInfo::GetEventUserInfo()) {
45 m_MDTHitColl = eventInfo->GetHitCollectionMap()->Find<MDTSimHitCollection>(m_hitCollectionName);
46 m_g4UserEventInfo = eventInfo;
47 }
49 // START OF COSMICS SPECIFIC CODE
50 m_mom = Amg::Vector3D(0.,0.,0.);
51 m_globH = Amg::Vector3D(0.,0.,0.);
52 // END OF COSMICS SPECIFIC CODE
53}
AtlasHitsVector< MDTSimHit > MDTSimHitCollection
static AtlasG4EventUserInfo * GetEventUserInfo()
Eigen::Matrix< double, 3, 1 > Vector3D

◆ ProcessHits()

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

Definition at line 55 of file MDTSensitiveDetectorCosmics.cxx.

55 {
56
57 if (!m_MDTHitColl) {
58 G4Exception("MDTSensitiveDetectorCosmics::ProcessHits", "MDTCosmicHitCollectionMissing", FatalException,
59 "Hit collection not initialized; did SetupEvent run?");
60 return false;
61 }
62
63 G4Track* currentTrack = aStep->GetTrack();
64
65 // MDTs sensitive to charged particle only
66 if (currentTrack->GetDefinition()->GetPDGCharge() == 0.0) {
67 if (currentTrack->GetDefinition()!=G4Geantino::GeantinoDefinition()) return true;
68 else if (currentTrack->GetDefinition()==G4ChargedGeantino::ChargedGeantinoDefinition()) return true;
69 }
70
71 G4VPhysicalVolume* physVolPostStep = aStep->GetPostStepPoint()->GetPhysicalVolume();
72 if (nullptr == physVolPostStep) return true;
73
74 // hit information to be recorded
75 double globalTime;
76 double driftRadius;
77 Amg::Vector3D localPosition;
78 localPosition.setZero();
79
80 // get top transformation
81 const G4AffineTransform trans = currentTrack->GetTouchable()->GetHistory()->GetTopTransform();
82
83 // START OF COSMICS SPECIFIC CODE
84 // global coordinates
85 G4ThreeVector globVrtx = aStep->GetPreStepPoint()->GetPosition();
86
87 // distance of the hit from (0,0,0) vertex - calculated from the PreStepPoint (approximation)
88 double dist = globVrtx.mag();
89 double lightspeed = 299.792458; /* in vacuo speed of light [mm/ns] */
90 double tOrigin = dist / lightspeed;
91
92 G4int trackid = aStep->GetTrack()->GetTrackID();
93 m_currVertex = Amg::Hep3VectorToEigen( aStep->GetTrack()->GetVertexPosition() );
94
95 // for cosmics: only primary muon tracks - track momentum when first entering the spectrometer (one muon per event)
96 if ((m_currVertex != m_vertex) && (trackid == 1)) {
97 // after calculationg the momentum magnidude, normalize it
98 m_mom = Amg::Hep3VectorToEigen( currentTrack->GetMomentum() );
99 m_momMag = m_mom.mag();
100 m_mom.normalize();
101 // the direction of the primary mu is used to calculate the t0, the position ot the t0, m_globH, is ONE for a track
102 Amg::Vector3D globVrtxFix = Amg::Hep3VectorToEigen( globVrtx );
103 double AlphaGlobal = -1*(globVrtxFix[0]*m_mom[0] + globVrtxFix[1]*m_mom[1] + globVrtxFix[2]*m_mom[2])/(m_mom[0]*m_mom[0] + m_mom[1]*m_mom[1] + m_mom[2]*m_mom[2]);
104 m_globH = globVrtxFix + AlphaGlobal*m_mom;
105// G4cout << "COSMICS MAIN TRACK IN THE MDT!" << G4endl;
106 }
107 double globalDist = sqrt((m_globH[0] - globVrtx[0])*(m_globH[0] - globVrtx[0]) +
108 (m_globH[1] - globVrtx[1])*(m_globH[1] - globVrtx[1]) +
109 (m_globH[2] - globVrtx[2])*(m_globH[2] - globVrtx[2]));
110 double tof = globalDist / lightspeed;
111 // END OF COSMICS SPECIFIC CODE
112
113 // transform pre and post step positions to local positions
114
115 Amg::Vector3D localVertex1( Amg::Hep3VectorToEigen( trans.TransformPoint(aStep->GetPreStepPoint()->GetPosition()) ) );
116 Amg::Vector3D localVertex2( Amg::Hep3VectorToEigen( trans.TransformPoint(aStep->GetPostStepPoint()->GetPosition()) ) );
117
118 // calculate local direction from begin- and end-point of the step
119 Amg::Vector3D localDirection( (localVertex2 - localVertex1) ); // normalized
120 localDirection.z() = 0.; // look in xy-plane
121
122 // See if particle passed wire by projecting begin- and end-point on the step direction
123 if( (localVertex2.dot(localDirection)) * (localVertex1.dot(localDirection)) < 0 ) { // particle passed wire
124
125 // compute drift radius ( = impact parameter)
126 double Xpos = localVertex1[0];
127 double Ypos = localVertex1[1];
128 double Xdir = localDirection[0];
129 double Ydir = localDirection[1];
130
131 double Alpha = -1*(Xpos*Xdir + Ypos*Ydir)/(Xdir*Xdir + Ydir*Ydir); // localPosition*localDirection
132
133 localPosition = localVertex1 + Alpha*(localVertex2-localVertex1);
134 driftRadius = localPosition.perp();
135 globalTime = aStep->GetPreStepPoint()->GetGlobalTime(); // take pre step time
136
137 }else{ // particle didn't pass wire, one of the end-points is the shortest distance to the wire
138
139 // calculate which of the two end-points is closer to the wire
140 double dist1 = localVertex1.perp();
141 double dist2 = localVertex2.perp();
142 if( dist1 < dist2 ){ // first is closer
143 driftRadius = dist1;
144 localPosition = localVertex1;
145 globalTime = aStep->GetPreStepPoint()->GetGlobalTime();
146 }else{ // second is closer
147 driftRadius = dist2;
148 localPosition = localVertex2;
149 globalTime = aStep->GetPostStepPoint()->GetGlobalTime();
150 }
151 }
152 m_vertex = Amg::Hep3VectorToEigen( aStep->GetTrack()->GetVertexPosition() );
153
154 if( driftRadius < m_driftRadius ){ // check if current step came closer to the wire than the previous ones
156 // if the track m_vertex is far from (0,0,0), takes the tof, otherwise take the "usual" g4 globalTime
157 ((((m_vertex.mag()) < 100) || ((fabs(globalTime - tOrigin)) < 0.1) ) ? (m_globalTime = globalTime) : (m_globalTime = tof));
158 // if m_globalTime != globalTime and m_globalTime != tof in the output, this is due to multiple hits
159 // before founding the good one (small approximation)
160 m_localPosition = localPosition;
161 }
162
163 // check if particle left tube or stopped in tube
164 G4String namePreStepMat = aStep->GetPreStepPoint()->GetMaterial()->GetName();
165 G4String namePostStepMat = aStep->GetPostStepPoint()->GetMaterial()->GetName();
166 G4String nameSD = aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetSensitiveDetector()->GetName();
167
168 // see if we were in the sensitive volume and left it, or the particle was killed
169 if( ((nameSD) && (namePreStepMat != namePostStepMat)) || (currentTrack->GetTrackStatus() == fStopAndKill)){
170
171 // get identifier
172 const G4TouchableHistory* touchHist = static_cast<const G4TouchableHistory*>(aStep->GetPreStepPoint()->GetTouchable());
173 int MDTid = GetIdentifier(touchHist);
174
175 TrackHelper trHelp(aStep->GetTrack());
176
177 // construct new mdt hit
179 trHelp.GenerateParticleLink(m_g4UserEventInfo ? m_g4UserEventInfo->GetEventStore() : nullptr),
180 aStep->GetStepLength(),
181 aStep->GetTotalEnergyDeposit(),
182 currentTrack->GetDefinition()->GetPDGEncoding(),
183 aStep->GetPreStepPoint()->GetKineticEnergy());
184
185 m_driftRadius = m_DEFAULT_TUBE_RADIUS; // reset start value of driftRadius
186 }
187 return true;
188}
int GetIdentifier(const G4TouchableHistory *touchHist)
Amg::Vector3D Hep3VectorToEigen(const CLHEP::Hep3Vector &CLHEPvector)
Converts a CLHEP-based CLHEP::Hep3Vector into an Eigen-based Amg::Vector3D.
@ driftRadius
trt, straws
Definition ParamDefs.h:53

Member Data Documentation

◆ m_currVertex

Amg::Vector3D MDTSensitiveDetectorCosmics::m_currVertex
private

Definition at line 108 of file MDTSensitiveDetectorCosmics.h.

◆ m_DEFAULT_TUBE_RADIUS

double MDTSensitiveDetectorCosmics::m_DEFAULT_TUBE_RADIUS
private

radius assigned to radius if radius is invalid

Definition at line 122 of file MDTSensitiveDetectorCosmics.h.

◆ m_driftRadius

double MDTSensitiveDetectorCosmics::m_driftRadius
private

Definition at line 117 of file MDTSensitiveDetectorCosmics.h.

◆ m_g4UserEventInfo

AtlasG4EventUserInfo* MDTSensitiveDetectorCosmics::m_g4UserEventInfo {nullptr}
private

Definition at line 114 of file MDTSensitiveDetectorCosmics.h.

114{nullptr};

◆ m_globalTime

double MDTSensitiveDetectorCosmics::m_globalTime
private

Definition at line 118 of file MDTSensitiveDetectorCosmics.h.

◆ m_globH

Amg::Vector3D MDTSensitiveDetectorCosmics::m_globH
private

Definition at line 109 of file MDTSensitiveDetectorCosmics.h.

◆ m_hitCollectionName

std::string MDTSensitiveDetectorCosmics::m_hitCollectionName
private

member data

Definition at line 112 of file MDTSensitiveDetectorCosmics.h.

◆ m_localPosition

Amg::Vector3D MDTSensitiveDetectorCosmics::m_localPosition
private

Definition at line 119 of file MDTSensitiveDetectorCosmics.h.

◆ m_MDTHitColl

MDTSimHitCollection* MDTSensitiveDetectorCosmics::m_MDTHitColl {nullptr}
private

Definition at line 113 of file MDTSensitiveDetectorCosmics.h.

113{nullptr};

◆ m_mom

Amg::Vector3D MDTSensitiveDetectorCosmics::m_mom
private

Definition at line 105 of file MDTSensitiveDetectorCosmics.h.

◆ m_momMag

double MDTSensitiveDetectorCosmics::m_momMag
private

Definition at line 106 of file MDTSensitiveDetectorCosmics.h.

◆ m_muonHelper

const MdtHitIdHelper* MDTSensitiveDetectorCosmics::m_muonHelper
private

Definition at line 115 of file MDTSensitiveDetectorCosmics.h.

◆ m_vertex

Amg::Vector3D MDTSensitiveDetectorCosmics::m_vertex
private

Definition at line 107 of file MDTSensitiveDetectorCosmics.h.


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