ATLAS Offline Software
Loading...
Searching...
No Matches
SctSensor_CTB.cxx
Go to the documentation of this file.
1/*
2 Copyright (C) 2002-2019 CERN for the benefit of the ATLAS collaboration
3*/
4
5//
6// SCT Sensitive Detector.
7// The Hits are processed here. For every hit I get the position and
8// an information on the sensor in which the interaction happened
9//
10//#define SCTSD_CTB_DEBUG
11
12// class header
13#include "SctSensor_CTB.h"
14
15// athena includes
17#include "MCTruth/TrackHelper.h"
18
19// Geant4 includes
20#include <G4EventManager.hh>
21#include "G4Step.hh"
22#include "G4ThreeVector.hh"
23#include "G4SDManager.hh"
24#include "G4Geantino.hh"
25#include "G4ChargedGeantino.hh"
26
27// CLHEP transform
28#include "CLHEP/Geometry/Transform3D.h"
29
30//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
31
32SctSensor_CTB::SctSensor_CTB(const std::string& name, const std::string& hitCollectionName)
33 : G4VSensitiveDetector(name)
34 , m_HitCollName(hitCollectionName)
35{
36}
37
38//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
39
40void SctSensor_CTB::Initialize(G4HCofThisEvent *)
41{
42 // ISF calls G4SDManager::PrepareNewEvent() before the Geant4 event loop starts...
43 if(auto* eventManger = G4EventManager::GetEventManager()){
44 if(auto* eventInfo = static_cast<AtlasG4EventUserInfo*>(eventManger->GetUserInformation())){
45 m_HitColl = eventInfo->GetHitCollectionMap()->Find<SiHitCollection>(m_HitCollName);
46 m_g4UserEventInfo = eventInfo;
47 }
48 }
49}
50
51//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
52
53G4bool SctSensor_CTB::ProcessHits(G4Step* aStep, G4TouchableHistory* /*ROhist*/)
54{
55 G4double edep = aStep->GetTotalEnergyDeposit();
56 if(edep==0.){
57 if(aStep->GetTrack()->GetDefinition()!=G4Geantino::GeantinoDefinition() &&
58 aStep->GetTrack()->GetDefinition()!=G4ChargedGeantino::ChargedGeantinoDefinition())
59 return false;
60 }
61 edep *= CLHEP::MeV;
62 //use the global time. i.e. the time from the beginning of the event
63 //
64 // Get the Touchable History:
65 //
66 const G4TouchableHistory* myTouch = dynamic_cast<const G4TouchableHistory*>(aStep->GetPreStepPoint()->GetTouchable());
67 if (not myTouch) {
68 G4cout << "SctSensor_CTB::ProcessHits bad dynamic_cast" << G4endl;
69 return false;
70 }
71 //
72 // Get the hit coordinates. Start and End Point
73 //
74 G4ThreeVector coord1 = aStep->GetPreStepPoint()->GetPosition();
75 G4ThreeVector coord2 = aStep->GetPostStepPoint()->GetPosition();
76 //
77 // Calculate the local step begin and end position.
78 // From a G4 FAQ:
79 // http://geant4-hn.slac.stanford.edu:5090/HyperNews/public/get/geometry/17/1.html
80 //
81 const G4AffineTransform transformation = myTouch->GetHistory()->GetTopTransform();
82 G4ThreeVector localPosition1 = transformation.TransformPoint(coord1);
83 G4ThreeVector localPosition2 = transformation.TransformPoint(coord2);
84
85 HepGeom::Point3D<double> lP1,lP2;
86 lP1[SiHit::xEta] = localPosition1[2]*CLHEP::mm;
87 lP1[SiHit::xPhi] = localPosition1[1]*CLHEP::mm;
88 lP1[SiHit::xDep] = localPosition1[0]*CLHEP::mm;
89
90 lP2[SiHit::xEta] = localPosition2[2]*CLHEP::mm;
91 lP2[SiHit::xPhi] = localPosition2[1]*CLHEP::mm;
92 lP2[SiHit::xDep] = localPosition2[0]*CLHEP::mm;
93
94 //
95 // Get it into a vector in local coords and with the right units:
96
97
98 // Now Navigate the history to know in what detector the step is:
99 // and finally set the ID of det element in which the hit is.
100 //
101 // G4int History;
102 //
103 // Is it TB, barrel or endcap?
104 //
105 int BrlEcap = 0;
106 int LayerDisk = 0;
107 int etaMod = 1;
108 int phiMod = 0;
109 int side = 0;
110 //
111 // In the case of the TB the positioning integers won't be initialized
112 // and the identifying integer will be zero. There is no need to do
113 // anything else
114
115 int BEcopyNo = myTouch->GetVolume()->GetCopyNo();
116
117#ifdef SCTSD_CTB_DEBUG
118 for (int i=0;i<myTouch->GetHistoryDepth();i++){
119 std::string detname=myTouch->GetVolume(i)->GetLogicalVolume()->GetName();
120 int copyno=myTouch->GetVolume(i)->GetCopyNo();
121 if (verboseLevel>1) G4cout << "Volume "<<detname <<" Copy Nr. "<<copyno << G4endl;
122 }
123#endif
124
125 // Testbeam geometry from GeoModel.
126 //
127 int BEcopyNoTest = BEcopyNo/100;
128 if(BEcopyNoTest == 5 || BEcopyNoTest == 6 || BEcopyNo == 1000 || BEcopyNo == 2000) {
129 BrlEcap=0;
130 etaMod=1;
131
132 if(BEcopyNoTest == 5 || BEcopyNoTest == 6 ) {
133 side = myTouch->GetVolume(0)->GetCopyNo()%2;
134 phiMod=myTouch->GetVolume(1)->GetCopyNo()%2;
135 LayerDisk = myTouch->GetVolume(1)->GetCopyNo()/2;
136 } else if ( BEcopyNo == 1000) {
137 side = myTouch->GetVolume(1)->GetCopyNo();
138 phiMod=myTouch->GetVolume(2)->GetCopyNo()%2;
139 LayerDisk = myTouch->GetVolume(2)->GetCopyNo()/2;
140 } // otherwise use side = phiMod = LayerDisk = 0.
141 // Some printout
142#ifdef SCTSD_CTB_DEBUG
143 if (verboseLevel>1) G4cout << "In the SCT TestBeam" << G4endl;
144 if (verboseLevel>1) G4cout << " * Module Number: " << LayerDisk << "," << phiMod << " Side: " << side << G4endl;
145 if (verboseLevel>5) G4cout << " Identifier will be [2.2."
146 << BrlEcap << "."
147 << LayerDisk << "."
148 << phiMod << "."
149 << etaMod << "."
150 << side << "]" << G4endl;
151#endif
152 } else {
153 G4ExceptionDescription description;
154 description << "ProcessHits: Unexpected copy number for sensor";
155 G4Exception("SctSensor_CTB", "UnexpectedCopyNumberForSctSensor", FatalException, description);
156 abort();
157 }
158 TrackHelper trHelp(aStep->GetTrack());
159 m_HitColl->Emplace(lP1,
160 lP2,
161 edep,
162 aStep->GetPreStepPoint()->GetGlobalTime(),
163 trHelp.GenerateParticleLink(m_g4UserEventInfo ? m_g4UserEventInfo->GetEventStore() : nullptr),
164 1,BrlEcap,LayerDisk,etaMod,phiMod,side);
165 return true;
166}
AtlasHitsVector< SiHit > SiHitCollection
This class is attached to G4Event objects as UserInformation.
G4bool ProcessHits(G4Step *, G4TouchableHistory *) override final
std::string m_HitCollName
SctSensor_CTB(const std::string &name, const std::string &hitCollectionName)
void Initialize(G4HCofThisEvent *) override final
AtlasG4EventUserInfo * m_g4UserEventInfo
SiHitCollection * m_HitColl
@ xPhi
Definition SiHit.h:162
@ xEta
Definition SiHit.h:162
@ xDep
Definition SiHit.h:162
HepMcParticleLink GenerateParticleLink()
Generates a creates new HepMcParticleLink object on the stack based on GetUniqueID(),...
Definition TrackHelper.h:52
std::string description
glabal timer - how long have I taken so far?
Definition hcg.cxx:91