ATLAS Offline Software
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 
32 SctSensor_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 
40 void 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  }
47  }
48 }
49 
50 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
51 
52 G4bool SctSensor_CTB::ProcessHits(G4Step* aStep, G4TouchableHistory* /*ROhist*/)
53 {
54  G4double edep = aStep->GetTotalEnergyDeposit();
55  if(edep==0.){
56  if(aStep->GetTrack()->GetDefinition()!=G4Geantino::GeantinoDefinition() &&
57  aStep->GetTrack()->GetDefinition()!=G4ChargedGeantino::ChargedGeantinoDefinition())
58  return false;
59  }
60  edep *= CLHEP::MeV;
61  //use the global time. i.e. the time from the beginning of the event
62  //
63  // Get the Touchable History:
64  //
65  const G4TouchableHistory* myTouch = dynamic_cast<const G4TouchableHistory*>(aStep->GetPreStepPoint()->GetTouchable());
66  if (not myTouch) {
67  G4cout << "SctSensor_CTB::ProcessHits bad dynamic_cast" << G4endl;
68  return false;
69  }
70  //
71  // Get the hit coordinates. Start and End Point
72  //
73  G4ThreeVector coord1 = aStep->GetPreStepPoint()->GetPosition();
74  G4ThreeVector coord2 = aStep->GetPostStepPoint()->GetPosition();
75  //
76  // Calculate the local step begin and end position.
77  // From a G4 FAQ:
78  // http://geant4-hn.slac.stanford.edu:5090/HyperNews/public/get/geometry/17/1.html
79  //
80  const G4AffineTransform transformation = myTouch->GetHistory()->GetTopTransform();
81  G4ThreeVector localPosition1 = transformation.TransformPoint(coord1);
82  G4ThreeVector localPosition2 = transformation.TransformPoint(coord2);
83 
84  HepGeom::Point3D<double> lP1,lP2;
85  lP1[SiHit::xEta] = localPosition1[2]*CLHEP::mm;
86  lP1[SiHit::xPhi] = localPosition1[1]*CLHEP::mm;
87  lP1[SiHit::xDep] = localPosition1[0]*CLHEP::mm;
88 
89  lP2[SiHit::xEta] = localPosition2[2]*CLHEP::mm;
90  lP2[SiHit::xPhi] = localPosition2[1]*CLHEP::mm;
91  lP2[SiHit::xDep] = localPosition2[0]*CLHEP::mm;
92 
93  //
94  // Get it into a vector in local coords and with the right units:
95 
96 
97  // Now Navigate the history to know in what detector the step is:
98  // and finally set the ID of det element in which the hit is.
99  //
100  // G4int History;
101  //
102  // Is it TB, barrel or endcap?
103  //
104  int BrlEcap = 0;
105  int LayerDisk = 0;
106  int etaMod = 1;
107  int phiMod = 0;
108  int side = 0;
109  //
110  // In the case of the TB the positioning integers won't be initialized
111  // and the identifying integer will be zero. There is no need to do
112  // anything else
113 
114  int BEcopyNo = myTouch->GetVolume()->GetCopyNo();
115 
116 #ifdef SCTSD_CTB_DEBUG
117  for (int i=0;i<myTouch->GetHistoryDepth();i++){
118  std::string detname=myTouch->GetVolume(i)->GetLogicalVolume()->GetName();
119  int copyno=myTouch->GetVolume(i)->GetCopyNo();
120  if (verboseLevel>1) G4cout << "Volume "<<detname <<" Copy Nr. "<<copyno << G4endl;
121  }
122 #endif
123 
124  // Testbeam geometry from GeoModel.
125  //
126  int BEcopyNoTest = BEcopyNo/100;
127  if(BEcopyNoTest == 5 || BEcopyNoTest == 6 || BEcopyNo == 1000 || BEcopyNo == 2000) {
128  BrlEcap=0;
129  etaMod=1;
130 
131  if(BEcopyNoTest == 5 || BEcopyNoTest == 6 ) {
132  side = myTouch->GetVolume(0)->GetCopyNo()%2;
133  phiMod=myTouch->GetVolume(1)->GetCopyNo()%2;
134  LayerDisk = myTouch->GetVolume(1)->GetCopyNo()/2;
135  } else if ( BEcopyNo == 1000) {
136  side = myTouch->GetVolume(1)->GetCopyNo();
137  phiMod=myTouch->GetVolume(2)->GetCopyNo()%2;
138  LayerDisk = myTouch->GetVolume(2)->GetCopyNo()/2;
139  } // otherwise use side = phiMod = LayerDisk = 0.
140  // Some printout
141 #ifdef SCTSD_CTB_DEBUG
142  if (verboseLevel>1) G4cout << "In the SCT TestBeam" << G4endl;
143  if (verboseLevel>1) G4cout << " * Module Number: " << LayerDisk << "," << phiMod << " Side: " << side << G4endl;
144  if (verboseLevel>5) G4cout << " Identifier will be [2.2."
145  << BrlEcap << "."
146  << LayerDisk << "."
147  << phiMod << "."
148  << etaMod << "."
149  << side << "]" << G4endl;
150 #endif
151  } else {
152  G4ExceptionDescription description;
153  description << "ProcessHits: Unexpected copy number for sensor";
154  G4Exception("SctSensor_CTB", "UnexpectedCopyNumberForSctSensor", FatalException, description);
155  abort();
156  }
157  TrackHelper trHelp(aStep->GetTrack());
158  m_HitColl->Emplace(lP1,
159  lP2,
160  edep,
161  aStep->GetPreStepPoint()->GetGlobalTime(),
162  trHelp.GenerateParticleLink(),
163  1,BrlEcap,LayerDisk,etaMod,phiMod,side);
164  return true;
165 }
AtlasG4EventUserInfo
This class is attached to G4Event objects as UserInformation. It holds a pointer to the HepMC::GenEve...
Definition: AtlasG4EventUserInfo.h:23
SctSensor_CTB::ProcessHits
G4bool ProcessHits(G4Step *, G4TouchableHistory *) override final
Definition: SctSensor_CTB.cxx:52
SiHit::xDep
@ xDep
Definition: SiHit.h:162
python.SystemOfUnits.mm
float mm
Definition: SystemOfUnits.py:98
SctSensor_CTB.h
TrackHelper.h
SiHit::xPhi
@ xPhi
Definition: SiHit.h:162
AtlasHitsVector< SiHit >
SctSensor_CTB::m_HitColl
SiHitCollection * m_HitColl
Definition: SctSensor_CTB.h:48
python.SystemOfUnits.MeV
float MeV
Definition: SystemOfUnits.py:172
TRT::Hit::side
@ side
Definition: HitInfo.h:83
AtlasHitsVector::Emplace
void Emplace(Args &&... args)
Definition: AtlasHitsVector.h:80
TrackHelper
Definition: TrackHelper.h:14
SctSensor_CTB::SctSensor_CTB
SctSensor_CTB(const std::string &name, const std::string &hitCollectionName)
Definition: SctSensor_CTB.cxx:32
LayerDisk
Definition: WaferTree.h:57
lumiFormat.i
int i
Definition: lumiFormat.py:85
SctSensor_CTB::Initialize
void Initialize(G4HCofThisEvent *) override final
Definition: SctSensor_CTB.cxx:40
name
std::string name
Definition: Control/AthContainers/Root/debug.cxx:240
SiHit::xEta
@ xEta
Definition: SiHit.h:162
AtlasG4EventUserInfo.h
TrackHelper::GenerateParticleLink
HepMcParticleLink GenerateParticleLink()
Generates a creates new HepMcParticleLink object on the stack based on GetUniqueID(),...
Definition: TrackHelper.h:35
fillSCTHists.etaMod
etaMod
Definition: fillSCTHists.py:23
description
std::string description
glabal timer - how long have I taken so far?
Definition: hcg.cxx:88
SctSensor_CTB::m_HitCollName
std::string m_HitCollName
Definition: SctSensor_CTB.h:47