ATLAS Offline Software
Loading...
Searching...
No Matches
RPCSensitiveDetector.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2024 CERN for the benefit of the ATLAS collaboration
3*/
4
6#include "G4ThreeVector.hh"
7#include "G4Trd.hh"
9#include <string>
10#include "G4Exception.hh"
11#include "G4Geantino.hh"
12#include "G4ChargedGeantino.hh"
13
14//#include "SimHelpers/DetectorGeometryHelper.h"
16#include "MCTruth/TrackHelper.h"
17#include <sstream>
18
19
22#include "GaudiKernel/SystemOfUnits.h"
23
24// construction/destruction
25RPCSensitiveDetector::RPCSensitiveDetector(const std::string& name, const std::string& hitCollectionName, unsigned int nGasGaps)
26 : G4VSensitiveDetector( name )
27 , m_hitCollectionName( hitCollectionName )
28{
30}
31
33{
34 m_myRPCHitColl = nullptr;
35 if (auto* eventInfo = AtlasG4EventUserInfo::GetEventUserInfo()) {
36 m_myRPCHitColl = eventInfo->GetHitCollectionMap()->Find<RPCSimHitCollection>(m_hitCollectionName);
37 m_g4UserEventInfo = eventInfo;
38 }
39 //FIXME probably only need to call this bit at start of the event
40 //loop rather than the start of each G4Event.
41 if (verboseLevel>5) G4cout << "Initializing SD" << G4endl;
42
43}
44
45G4bool RPCSensitiveDetector::ProcessHits(G4Step* aStep,G4TouchableHistory*) {
46
47 if (!m_myRPCHitColl) {
48 G4Exception("RPCSensitiveDetector::ProcessHits", "RPCHitCollectionMissing", FatalException,
49 "Hit collection not initialized; did SetupEvent run?");
50 return false;
51 }
52
53 G4Track* track = aStep->GetTrack();
54
56 if (track->GetDefinition()->GetPDGCharge() == 0.0) {
57 if (track->GetDefinition()!=G4Geantino::GeantinoDefinition()) {
58 return true;
59 }
60 }
61
62 const G4TouchableHistory* touchHist = static_cast<const G4TouchableHistory*>(aStep->GetPreStepPoint()->GetTouchable());
63 G4ThreeVector position = aStep->GetPreStepPoint()->GetPosition();
64 G4ThreeVector postPosition = aStep->GetPostStepPoint()->GetPosition();
65 const G4AffineTransform trans = track->GetTouchable()->GetHistory()->GetTopTransform(); // from global to local
66
67 // necessary to assign correct identifiers
68 int rpcIsRotated = 0;
69 // int layer = 0;
70
71 std::string tech;
72 // fields for the RPC identifier construction
73 std::string stationName;
74 int stationEta= 0;
75 int stationPhi= 0;
76 int doubletR= 0;
77 // int doubletZ= 0;
78 int doubletPhi= 0;
79 int gasGap= 0;
80 // int eta_inversed=0;
81 // int phi_inversed=0;
82 int station_rotated=0; // tells us if the station was rotated before being positioned.
83 // if =1 we have to correct some IDs
84 int zNeg_original=0; // tells if the station at z<0 was obtained duplicating a station at z>0
85
86 // RPC hit information
87 double globalTime = aStep->GetPreStepPoint()->GetGlobalTime();
88 Amg::Vector3D localPosition = Amg::Hep3VectorToEigen( trans.TransformPoint(position) );
89 Amg::Vector3D localPostPosition = Amg::Hep3VectorToEigen( trans.TransformPoint(postPosition) );
90 {
91 const Amg::Vector3D stepVector = localPostPosition - localPosition;
93 if (stepVector.mag()>std::numeric_limits<float>::epsilon() &&
94 std::abs(std::abs(Amg::angle(stepVector, Amg::Vector3D::UnitX()))
95 - 90.*Gaudi::Units::deg) < 0.0001* Gaudi::Units::deg) {
96 return true;
97 }
98 }
99 int mydbZ=0;
100 int mydbPMod=0;
101 int mydbP=0;
102 bool isAssembly = false;
103 // scan geometry tree to identify the hit channel
104 for (int i=touchHist->GetHistoryDepth();i>=0;i--) {
105
106 std::string::size_type npos;
107 std::string volName = touchHist->GetVolume(i)->GetName();
108 std::string num=volName.substr(3,2);
109 if(num[0]==' ') num[0]=0;
110
111 // check if this station is an assembly
112 if ((npos = volName.find("av_")) != std::string::npos &&
113 (npos = volName.find("impr_")) != std::string::npos) isAssembly = true;
114
115 // stationName, stationEta, stationPhi
116 if ((npos = volName.find("station")) != std::string::npos && (!isAssembly)) {
117
118 stationName = volName.substr(0,npos-1);
119
120 int volCopyNo = touchHist->GetVolume(i)->GetCopyNo();
121
122 if(abs(volCopyNo/1000)==1){
123 zNeg_original=1;
124 volCopyNo=volCopyNo%1000;
125 }
126
127 stationEta = volCopyNo/100;
128 stationPhi = abs(volCopyNo%100);
129
130 if(stationEta<0&&!zNeg_original) station_rotated=1;
131
132 // stationName, stationEta, stationPhi
133 } else if ((npos = volName.find("RPC")) != std::string::npos && isAssembly) {
134 // vol name for Assembly components are
135 // av_WWW_impr_XXX_Muon::BMSxMDTxx_pv_ZZZ_NAME
136 // where WWW is ass. istance nr.
137 // XXX is comp. imprint nr.
138 // BMSxMDTxx is the name of the comp. log.Vol.
139 // x station sub-type; xx technology subtype
140 // ZZZ is the comp. number of order
141 // NAME is the comp. tag (geoIdentifierTag)
142 // for RPCs it is SwapdbPdbRdbZ[aaa]
143 // with aaa = doubletZ + doubletR*100 + dbphi*1000;
144 // aaa = -aaa if iswap == -1
145 // dbphi = 1 always but for S1 or S2 at doubletPhi = 2
146 // copy numbers for Ass.components are =
147 // CopyNoBase(= geoIdentifierTag of the assembly) + geoIdentifierTag of the component
148 // geoIdentifierTag of the component = Job
149 // geoIdentifierTag of the assembly = (sideC*10000 +
150 // mirsign*1000 + abs(zi)*100 + fi+1)*100000;
151 // mirsign = 1 if zi<0 and !mirrored; otherwise 0
152 std::string::size_type loc1,loc2;
153 if ((loc1 = volName.find("Muon::")) != std::string::npos) {
154 stationName = volName.substr(loc1+6,4); //type and subtype
155 }
156
157 int volCopyNo = touchHist->GetVolume(i)->GetCopyNo();
158 int copyNrBase = int(volCopyNo/100000);
159 int sideC = int(copyNrBase/10000);
160 int zi = int((copyNrBase%1000)/100);
161 int mirfl = int((copyNrBase%10000)/1000);
162 int fi = int(copyNrBase%100);
163 if (sideC == 1) zi = -zi;
164 stationEta = zi;
165 stationPhi = fi;
166 zNeg_original = mirfl;
167
168 if(stationEta<0&&!zNeg_original) station_rotated=1;
169
170 // doubletZ = 1;
171 tech=volName.substr(npos,5);
172
173 // now get the geoIdentifierTag of the rpc components
174 int gmID = 0;
175 if ((loc1 = volName.find('[')) != std::string::npos) {
176 if ((loc2 = volName.find(']', loc1+1)) != std::string::npos) {
177 //G4cout << "first [ is at "<<loc1<<" first ] at "<<loc2 << G4endl;
178 std::istringstream istrvar(volName.substr(loc1+1,loc2-loc1-1));
179 istrvar>>gmID;
180 }
181 }
182 int kk=gmID;
183
184 if (kk < 0) rpcIsRotated=1;
185
186 doubletR =(abs(kk)%1000)/100;
187 mydbZ = abs(int(kk%10));
188 mydbPMod = abs(int(kk/1000));
189 // doubleR, doubletZ
190 } else if ((npos = volName.find("rpccomponent")) != std::string::npos && (!isAssembly)) {
191
192 std::string::size_type loc1,loc2;
193 tech=volName.substr(npos-5,5);
194 int gmID = 0;
195 if ((loc1 = volName.find('[')) != std::string::npos) {
196 if ((loc2 = volName.find(']', loc1+1)) != std::string::npos) {
197 std::istringstream istrvar(volName.substr(loc1+1,loc2-loc1-1));
198 istrvar>>gmID;
199 }
200 }
201 mydbZ = abs(int(gmID%10));
202 mydbPMod = abs(int(gmID/1000));
203
204 int kk=touchHist->GetVolume(i)->GetCopyNo();
205
206 if (kk < 0) rpcIsRotated=1;
207
208 doubletR=(abs(kk)%1000)/100;
209
210 } else if ((npos = volName.find("layer")) != std::string::npos) {
211
212 int copyNo = touchHist->GetVolume(i)->GetCopyNo();
213
214 if (copyNo == 1) {
215 rpcIsRotated ? gasGap = 2 : gasGap = 1;
216 } else if (copyNo ==2) {
217 rpcIsRotated ? gasGap = 1 : gasGap = 2;
218 } else if (copyNo ==3) {
219 gasGap = 3;
220 }
221 } else if((npos = volName.find("gas volume")) != std::string::npos) {
222
223 int copyNo = touchHist->GetVolume(i)->GetCopyNo();
224 if (copyNo == 0) {
225 doubletPhi = 1;
226 } else if (copyNo == 10) {
227 doubletPhi = 2;
228 }
229 mydbP = doubletPhi;
230 int ngap_in_s=0;
231 int nstrippanel_in_s=0;
232 std::string::size_type loc1;
233 if ((loc1 = volName.find("gg_in_s")) != std::string::npos) {
234 std::istringstream istrvar(volName.substr(loc1-1,1));
235 istrvar>>ngap_in_s;
236 }
237 if ((loc1 = volName.find("sp_in_s")) != std::string::npos) {
238 std::istringstream istrvar(volName.substr(loc1-1,1));
239 istrvar>>nstrippanel_in_s;
240 }
241 if (ngap_in_s == 1 && nstrippanel_in_s == 2) {
242 if(localPosition.y()>0) mydbP=2;
243 else mydbP=1;
244 } else if (ngap_in_s == 1 && nstrippanel_in_s == 1) {
245 mydbP = mydbPMod;
246 }
247 }
248 }
249
251 // correct wrong IDs due to station rotation. only doubletZ and doubletPhi are affected at this point
253
254 if(station_rotated){
255 // doubletPhi correction
256 // note that this is correct also for ribs chambers
257 mydbP++;
258 if (mydbP>2) mydbP=1;
259 // G4cout << "Rcd PC - swap dbPhi because ch is MIRRORED "<<doubletPhi << G4endl;
260
261 // strip numbering: if the station is rotated both eta and phi directions get inversed.
262 // commented out for geomodel!
263 }
264
266
267 // now we have the position in the gas gap, with correct axis orientation, and the offlineIDs *of the strip panel*
268 // nstrip is supposed to be calculated by RPC_Digitizer.
269
270 // construct two (one in eta and one in phi) new RPC hits and store them in hit collection
271
272 //construct the hit identifiers
273 if (verboseLevel>5) {
274 G4cout << "hit in station "<<stationName<< " on technology "<<tech << G4endl;
275 G4cout << "constructing ids (stName, stEta, stPhi, dr, dZ, dPhi)= "<<stationName<< " "<< stationEta<<" " << stationPhi<< " "<<doubletR<< " "<< mydbZ<< " "<<mydbP << G4endl;
276 }
277
278 HitID RPCid_eta = m_muonHelper->BuildRpcHitId(stationName, stationPhi, stationEta,
279 mydbZ, doubletR, gasGap, mydbP,0);
280
281 HitID RPCid_phi = m_muonHelper->BuildRpcHitId(stationName, stationPhi, stationEta,
282 mydbZ, doubletR, gasGap, mydbP,1);
283
284 // retrieve track barcode
285 TrackHelper trHelp(aStep->GetTrack());
286
287 //construct new rpc hit
288 m_myRPCHitColl->Emplace(RPCid_eta, globalTime,
289 localPosition,
290 trHelp.GenerateParticleLink(m_g4UserEventInfo ? m_g4UserEventInfo->GetEventStore() : nullptr),
291 localPostPosition,
292 aStep->GetTotalEnergyDeposit(),
293 aStep->GetStepLength(),
294 track->GetDefinition()->GetPDGEncoding(),
295 aStep->GetPreStepPoint()->GetKineticEnergy());
296 m_myRPCHitColl->Emplace(RPCid_phi, globalTime,
297 localPosition,
298 trHelp.GenerateParticleLink(m_g4UserEventInfo ? m_g4UserEventInfo->GetEventStore() : nullptr),
299 localPostPosition,
300 aStep->GetTotalEnergyDeposit(),
301 aStep->GetStepLength(),
302 track->GetDefinition()->GetPDGEncoding(),
303 aStep->GetPreStepPoint()->GetKineticEnergy());
304
305 return true;
306}
int HitID
AtlasHitsVector< RPCSimHit > RPCSimHitCollection
static AtlasG4EventUserInfo * GetEventUserInfo()
RPCSensitiveDetector(const std::string &name, const std::string &hitCollectionName, unsigned int nGasGaps)
construction/destruction
std::string m_hitCollectionName
member data
G4bool ProcessHits(G4Step *, G4TouchableHistory *) override final
const RpcHitIdHelper * m_muonHelper
void Initialize(G4HCofThisEvent *) override final
member functions
AtlasG4EventUserInfo * m_g4UserEventInfo
RPCSimHitCollection * m_myRPCHitColl
static const RpcHitIdHelper * GetHelper(unsigned int nGasGaps=2)
HepMcParticleLink GenerateParticleLink()
Generates a creates new HepMcParticleLink object on the stack based on GetUniqueID(),...
Definition TrackHelper.h:52
double angle(const Amg::Vector3D &v1, const Amg::Vector3D &v2)
calculates the opening angle between two vectors
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