CMS 3D CMS Logo

SensitiveDetector.cc
Go to the documentation of this file.
2 
6 
7 #include "G4SDManager.hh"
8 #include "G4Step.hh"
9 #include "G4Track.hh"
10 #include "G4StepPoint.hh"
11 #include "G4Transform3D.hh"
12 #include "G4LogicalVolumeStore.hh"
13 #include "G4TouchableHistory.hh"
14 #include "G4VUserTrackInformation.hh"
15 
16 #include <sstream>
17 
19  const DDCompactView & cpv,
20  const SensitiveDetectorCatalog & clg,
21  edm::ParameterSet const & p) :
22  G4VSensitiveDetector(iname)
23 {
24  // for CMS hits
25  namesOfSD.push_back(iname);
26 
27  // Geant4 hit collection
28  collectionName.insert(iname);
29 
30  // register sensitive detector
31  G4SDManager * SDman = G4SDManager::GetSDMpointer();
32  SDman->AddNewDetector(this);
33 
34  const std::vector<std::string>& lvNames = clg.logicalNames(iname);
35  std::stringstream ss;
36  for (auto & lvname : lvNames) {
37  this->AssignSD(lvname);
38  ss << " " << lvname;
39  }
40  edm::LogInfo("SensitiveDetector") << " <" << iname <<"> : Assigns SD to LVs "
41  << ss.str();
42 
43 }
44 
46 
47 void SensitiveDetector::Initialize(G4HCofThisEvent * eventHC) {}
48 
49 void SensitiveDetector::EndOfEvent(G4HCofThisEvent * eventHC) {}
50 
52 {
53  G4LogicalVolumeStore * theStore = G4LogicalVolumeStore::GetInstance();
54  for (auto & lv : *theStore)
55  {
56  if (vname == lv->GetName()) {
57  lv->SetSensitiveDetector(this);
58  }
59  }
60 }
61 
63 {
64  const G4StepPoint * preStepPoint = step->GetPreStepPoint();
65  const G4ThreeVector& globalCoordinates = preStepPoint->GetPosition();
66  if (cd == WorldCoordinates) { return ConvertToLocal3DPoint(globalCoordinates); }
67  G4TouchableHistory * theTouchable=(G4TouchableHistory *)(preStepPoint->GetTouchable());
68  const G4ThreeVector localCoordinates = theTouchable->GetHistory()
69  ->GetTopTransform().TransformPoint(globalCoordinates);
70  return ConvertToLocal3DPoint(localCoordinates);
71 }
72 
74 {
75  const G4StepPoint * postStepPoint = step->GetPostStepPoint();
76  const G4ThreeVector& globalCoordinates = postStepPoint->GetPosition();
77  if (cd == WorldCoordinates) { return ConvertToLocal3DPoint(globalCoordinates); }
78  const G4StepPoint * preStepPoint = step->GetPreStepPoint();
79  G4TouchableHistory * theTouchable = (G4TouchableHistory *)(preStepPoint->GetTouchable());
80  const G4ThreeVector localCoordinates = theTouchable->GetHistory()
81  ->GetTopTransform().TransformPoint(globalCoordinates);
82  return ConvertToLocal3DPoint(localCoordinates);
83 }
84 
86 {
87  const G4StepPoint * preStepPoint = step->GetPreStepPoint();
88  G4TouchableHistory * theTouchable=(G4TouchableHistory *)(preStepPoint->GetTouchable());
89  G4ThreeVector localCoordinates = theTouchable->GetHistory()
90  ->GetTopTransform().TransformPoint(preStepPoint->GetPosition());
91  return ConvertToLocal3DPoint(localCoordinates);
92 }
93 
95 {
96  const G4ThreeVector& globalCoordinates = step->GetPostStepPoint()->GetPosition();
97  G4TouchableHistory * theTouchable =
98  (G4TouchableHistory *)(step->GetPreStepPoint()->GetTouchable());
99  G4ThreeVector localCoordinates = theTouchable->GetHistory()
100  ->GetTopTransform().TransformPoint(globalCoordinates);
101  return ConvertToLocal3DPoint(localCoordinates);
102 }
103 
105 {
106  TrackInformation* info = (TrackInformation*)(aTrack->GetUserInformation());
107  if (!info) {
108  edm::LogWarning("SensitiveDetector")
109  << " no TrackInformation available for trackID= "
110  << aTrack->GetTrackID();
111  throw SimG4Exception("SimG4CoreSensitiveDetector: cannot handle hits for "
112  + GetName());
113  }
114  return info;
115 }
116 
117 void SensitiveDetector::setNames(const std::vector<std::string>& hnames)
118 {
119  namesOfSD.clear();
120  namesOfSD = hnames;
121 }
122 
123 void SensitiveDetector::NaNTrap(const G4Step* aStep) const
124 {
125  const G4Track* CurrentTrk = aStep->GetTrack();
126  const G4ThreeVector& CurrentPos = CurrentTrk->GetPosition();
127  double xyz = CurrentPos.x() + CurrentPos.y() + CurrentPos.z();
128  const G4ThreeVector& CurrentMom = CurrentTrk->GetMomentum();
129  xyz += CurrentMom.x() + CurrentMom.y() + CurrentMom.z();
130 
131  if( edm::isNotFinite(xyz) ) {
132 
133  const G4VPhysicalVolume* pCurrentVol = aStep->GetPreStepPoint()->GetPhysicalVolume();
134  const G4String& NameOfVol = pCurrentVol->GetName();
135  throw SimG4Exception("SimG4CoreSensitiveDetector: Corrupted Event - NaN detected in volume "
136  + NameOfVol);
137  }
138 }
static const TGPicture * info(bool iBackgroundIsBlack)
SensitiveDetector(const std::string &iname, const DDCompactView &cpv, const SensitiveDetectorCatalog &, edm::ParameterSet const &p)
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
~SensitiveDetector() override
const std::vector< std::string > & logicalNames(const std::string &readoutName) const
Compact representation of the geometrical detector hierarchy.
Definition: DDCompactView.h:83
void AssignSD(const std::string &vname)
bool isNotFinite(T x)
Definition: isFinite.h:10
std::string const collectionName[nCollections]
Definition: Collections.h:47
Local3DPoint FinalStepPosition(const G4Step *step, coordinates) const
TrackInformation * cmsTrackInformation(const G4Track *aTrack)
std::vector< std::string > namesOfSD
void EndOfEvent(G4HCofThisEvent *eventHC) override
Local3DPoint LocalPostStepPosition(const G4Step *step) const
void setNames(const std::vector< std::string > &)
Local3DPoint LocalPreStepPosition(const G4Step *step) const
step
void Initialize(G4HCofThisEvent *eventHC) override
void NaNTrap(const G4Step *step) const
Local3DPoint InitialStepPosition(const G4Step *step, coordinates) const