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 edm::EventSetup& es,
20  const SensitiveDetectorCatalog& clg,
21  edm::ParameterSet const& p,
22  bool calo)
23  : G4VSensitiveDetector(iname), m_isCalo(calo) {
24  // for CMS hits
25  m_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_view>& lvNames = clg.logicalNames(iname);
35  std::stringstream ss;
36  for (auto& lvname : lvNames) {
37  this->AssignSD({lvname.data(), lvname.size()});
38  ss << " " << lvname;
39  }
40  edm::LogVerbatim("SensitiveDetector") << " <" << iname << "> : Assigns SD to LVs " << ss.str();
41 }
42 
44 
45 void SensitiveDetector::Initialize(G4HCofThisEvent* eventHC) {}
46 
47 void SensitiveDetector::EndOfEvent(G4HCofThisEvent* eventHC) {}
48 
50  G4LogicalVolumeStore* theStore = G4LogicalVolumeStore::GetInstance();
51  for (auto& lv : *theStore) {
52  if (vname == lv->GetName()) {
53  lv->SetSensitiveDetector(this);
54  }
55  }
56 }
57 
59  const G4StepPoint* preStepPoint = step->GetPreStepPoint();
60  const G4ThreeVector& globalCoordinates = preStepPoint->GetPosition();
61  if (cd == WorldCoordinates) {
62  return ConvertToLocal3DPoint(globalCoordinates);
63  }
64  const G4TouchableHistory* theTouchable = static_cast<const G4TouchableHistory*>(preStepPoint->GetTouchable());
65  const G4ThreeVector localCoordinates =
66  theTouchable->GetHistory()->GetTopTransform().TransformPoint(globalCoordinates);
67  return ConvertToLocal3DPoint(localCoordinates);
68 }
69 
71  const G4StepPoint* postStepPoint = step->GetPostStepPoint();
72  const G4ThreeVector& globalCoordinates = postStepPoint->GetPosition();
73  if (cd == WorldCoordinates) {
74  return ConvertToLocal3DPoint(globalCoordinates);
75  }
76  const G4StepPoint* preStepPoint = step->GetPreStepPoint();
77  const G4ThreeVector localCoordinates =
78  preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(globalCoordinates);
79  return ConvertToLocal3DPoint(localCoordinates);
80 }
81 
83  const G4StepPoint* preStepPoint = step->GetPreStepPoint();
84  G4ThreeVector localCoordinates =
85  preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(preStepPoint->GetPosition());
86  return ConvertToLocal3DPoint(localCoordinates);
87 }
88 
90  const G4ThreeVector& globalCoordinates = step->GetPostStepPoint()->GetPosition();
91  G4ThreeVector localCoordinates =
92  step->GetPreStepPoint()->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(globalCoordinates);
93  return ConvertToLocal3DPoint(localCoordinates);
94 }
95 
97  TrackInformation* info = (TrackInformation*)(aTrack->GetUserInformation());
98  if (!info) {
99  edm::LogWarning("SensitiveDetector") << " no TrackInformation available for trackID= " << aTrack->GetTrackID();
100  throw SimG4Exception("SimG4CoreSensitiveDetector: cannot handle hits for " + GetName());
101  }
102  return info;
103 }
104 
105 void SensitiveDetector::setNames(const std::vector<std::string>& hnames) {
106  m_namesOfSD.clear();
107  m_namesOfSD = hnames;
108 }
109 
110 void SensitiveDetector::NaNTrap(const G4Step* aStep) const {
111  const G4Track* CurrentTrk = aStep->GetTrack();
112  const G4ThreeVector& CurrentPos = CurrentTrk->GetPosition();
113  double xyz = CurrentPos.x() + CurrentPos.y() + CurrentPos.z();
114  const G4ThreeVector& CurrentMom = CurrentTrk->GetMomentum();
115  xyz += CurrentMom.x() + CurrentMom.y() + CurrentMom.z();
116 
117  if (edm::isNotFinite(xyz)) {
118  const G4VPhysicalVolume* pCurrentVol = aStep->GetPreStepPoint()->GetPhysicalVolume();
119  const G4String& NameOfVol = pCurrentVol->GetName();
120  throw SimG4Exception("SimG4CoreSensitiveDetector: Corrupted Event - NaN detected in volume " + NameOfVol);
121  }
122 }
static const TGPicture * info(bool iBackgroundIsBlack)
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
~SensitiveDetector() override
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
const std::vector< std::string_view > logicalNames(const std::string &readoutName) const
void AssignSD(const std::string &vname)
std::vector< std::string > m_namesOfSD
Local3DPoint FinalStepPosition(const G4Step *step, coordinates) const
TrackInformation * cmsTrackInformation(const G4Track *aTrack)
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
Definition: StallMonitor.cc:94
SensitiveDetector(const std::string &iname, const edm::EventSetup &es, const SensitiveDetectorCatalog &, edm::ParameterSet const &p, bool calo)
void Initialize(G4HCofThisEvent *eventHC) override
void NaNTrap(const G4Step *step) const
Local3DPoint InitialStepPosition(const G4Step *step, coordinates) const