ATLAS Offline Software
Public Member Functions | Private Member Functions | Private Attributes | List of all members
PixelSensorGmxSD Class Reference

#include <PixelSensorGmxSD.h>

Inheritance diagram for PixelSensorGmxSD:
Collaboration diagram for PixelSensorGmxSD:

Public Member Functions

 PixelSensorGmxSD (const std::string &name, const std::string &hitCollectionName, GeoModelIO::ReadGeoModel *sqlreader=nullptr)
 
virtual ~PixelSensorGmxSD ()
 
virtual G4bool ProcessHits (G4Step *, G4TouchableHistory *) override final
 
virtual void Initialize (G4HCofThisEvent *) override final
 
template<class... Args>
void AddHit (Args &&... args)
 Templated method to stuff a single hit into the sensitive detector class. More...
 

Private Member Functions

 FRIEND_TEST (PixelSensorGmxSDtest, Initialize)
 
 FRIEND_TEST (PixelSensorGmxSDtest, ProcessHits)
 
 FRIEND_TEST (PixelSensorGmxSDtest, AddHit)
 

Private Attributes

SG::WriteHandle< SiHitCollectionm_HitColl
 
GeoModelIO::ReadGeoModel * m_sqlreader
 

Detailed Description

Definition at line 28 of file PixelSensorGmxSD.h.

Constructor & Destructor Documentation

◆ PixelSensorGmxSD()

PixelSensorGmxSD::PixelSensorGmxSD ( const std::string &  name,
const std::string &  hitCollectionName,
GeoModelIO::ReadGeoModel *  sqlreader = nullptr 
)

Definition at line 37 of file PixelSensorGmxSD.cxx.

38  : G4VSensitiveDetector( name )
39  , m_HitColl( hitCollectionName )
40 {
41  m_sqlreader = sqlreader;
42 }

◆ ~PixelSensorGmxSD()

virtual PixelSensorGmxSD::~PixelSensorGmxSD ( )
inlinevirtual

Definition at line 38 of file PixelSensorGmxSD.h.

38 { /* I don't own myHitColl if all has gone well */ }

Member Function Documentation

◆ AddHit()

template<class... Args>
void PixelSensorGmxSD::AddHit ( Args &&...  args)
inline

Templated method to stuff a single hit into the sensitive detector class.

This could get rather tricky, but the idea is to allow fast simulations to use the very same SD classes as the standard simulation.

Definition at line 49 of file PixelSensorGmxSD.h.

49 { m_HitColl->Emplace( args... ); }

◆ FRIEND_TEST() [1/3]

PixelSensorGmxSD::FRIEND_TEST ( PixelSensorGmxSDtest  ,
AddHit   
)
private

◆ FRIEND_TEST() [2/3]

PixelSensorGmxSD::FRIEND_TEST ( PixelSensorGmxSDtest  ,
Initialize   
)
private

◆ FRIEND_TEST() [3/3]

PixelSensorGmxSD::FRIEND_TEST ( PixelSensorGmxSDtest  ,
ProcessHits   
)
private

◆ Initialize()

void PixelSensorGmxSD::Initialize ( G4HCofThisEvent *  )
finaloverridevirtual

Definition at line 45 of file PixelSensorGmxSD.cxx.

46 {
47  if (!m_HitColl.isValid()) m_HitColl = std::make_unique<SiHitCollection>();
48 }

◆ ProcessHits()

G4bool PixelSensorGmxSD::ProcessHits ( G4Step *  aStep,
G4TouchableHistory *   
)
finaloverridevirtual

Definition at line 51 of file PixelSensorGmxSD.cxx.

52 {
53  if (verboseLevel>5) G4cout << "Process Hit" << G4endl;
54 
55  G4double edep = aStep->GetTotalEnergyDeposit();
56  edep *= CLHEP::MeV;
57  if(edep==0.) {
58  if(aStep->GetTrack()->GetDefinition() != G4Geantino::GeantinoDefinition() &&
59  aStep->GetTrack()->GetDefinition() != G4ChargedGeantino::ChargedGeantinoDefinition())
60  return false;
61  }
62 
63  //use the global time. i.e. the time from the beginning of the event
64  //
65  // Get the Touchable History:
66  //
67  const G4TouchableHistory *myTouch = dynamic_cast<const G4TouchableHistory*>(aStep->GetPreStepPoint()->GetTouchable());
68  if (not myTouch) {
69  G4cout << "PixelSensorGmxSD::ProcessHits bad dynamic_cast" << G4endl;
70  return false;
71  }
72  if(verboseLevel>5){
73  for (int i=0;i<myTouch->GetHistoryDepth();i++){
74  std::string detname=myTouch->GetVolume(i)->GetLogicalVolume()->GetName();
75  int copyno=myTouch->GetVolume(i)->GetCopyNo();
76  G4cout << "Volume " <<detname <<" Copy Nr. " << copyno << G4endl;
77  }
78  }
79  //
80  // Get the hit coordinates. Start and End Point
81  //
82  G4ThreeVector coord1 = aStep->GetPreStepPoint()->GetPosition();
83  G4ThreeVector coord2 = aStep->GetPostStepPoint()->GetPosition();
84 
85  // Calculate the local step begin and end position.
86  // From a G4 FAQ:
87  // http://geant4-hn.slac.stanford.edu:5090/HyperNews/public/get/geometry/17/1.html
88  //
89  const G4AffineTransform transformation = myTouch->GetHistory()->GetTopTransform();
90  G4ThreeVector localPosition1 = transformation.TransformPoint(coord1);
91  G4ThreeVector localPosition2 = transformation.TransformPoint(coord2);
92 
93  HepGeom::Point3D<double> lP1,lP2;
94  lP1[SiHit::xEta] = localPosition1[2]*CLHEP::mm;
95  lP1[SiHit::xPhi] = localPosition1[1]*CLHEP::mm;
96  lP1[SiHit::xDep] = localPosition1[0]*CLHEP::mm;
97 
98  lP2[SiHit::xEta] = localPosition2[2]*CLHEP::mm;
99  lP2[SiHit::xPhi] = localPosition2[1]*CLHEP::mm;
100  lP2[SiHit::xDep] = localPosition2[0]*CLHEP::mm;
101 
102  TrackHelper trHelp(aStep->GetTrack());
103 
104  if(m_sqlreader){
105  //if sqlite inputs, Identifier indices come from PhysVol Name
106  std::string physVolName = myTouch->GetVolume(0)->GetName();
107 
108  int hitIdOfWafer = SiHitIdHelper::GetHelper()->buildHitIdFromStringITk(0,physVolName);
109 
110 
111  m_HitColl->Emplace(lP1,
112  lP2,
113  edep,
114  aStep->GetPreStepPoint()->GetGlobalTime(),//use the global time. i.e. the time from the beginning of the event
115  trHelp.GenerateParticleLink(),
116  hitIdOfWafer);
117  return true;
118 
119  }
120  // if not from SQLite, we assume that the Identifier has already been written in as the copy number
121  // (it should have done if GeoModel building ran within Athena)
122 
123  int id = myTouch->GetVolume()->GetCopyNo();
124 
125  m_HitColl->Emplace(lP1,
126  lP2,
127  edep,
128  aStep->GetPreStepPoint()->GetGlobalTime(),
129  trHelp.GenerateParticleLink(),
130  id);
131  return true;
132 }

Member Data Documentation

◆ m_HitColl

SG::WriteHandle<SiHitCollection> PixelSensorGmxSD::m_HitColl
private

Definition at line 53 of file PixelSensorGmxSD.h.

◆ m_sqlreader

GeoModelIO::ReadGeoModel* PixelSensorGmxSD::m_sqlreader
private

Definition at line 54 of file PixelSensorGmxSD.h.


The documentation for this class was generated from the following files:
SiHit::xPhi
@ xPhi
Definition: SiHit.h:162
python.SystemOfUnits.MeV
int MeV
Definition: SystemOfUnits.py:154
TrackHelper
Definition: TrackHelper.h:14
lumiFormat.i
int i
Definition: lumiFormat.py:92
SiHitIdHelper::buildHitIdFromStringITk
int buildHitIdFromStringITk(int part, const std::string &) const
Definition: SiHitIdHelper.cxx:131
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
SiHitIdHelper::GetHelper
static const SiHitIdHelper * GetHelper()
Definition: SiHitIdHelper.cxx:19
PixelSensorGmxSD::m_HitColl
SG::WriteHandle< SiHitCollection > m_HitColl
Definition: PixelSensorGmxSD.h:53
PixelSensorGmxSD::m_sqlreader
GeoModelIO::ReadGeoModel * m_sqlreader
Definition: PixelSensorGmxSD.h:54
python.CaloScaleNoiseConfig.args
args
Definition: CaloScaleNoiseConfig.py:80