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