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  bool calo) :
23  G4VSensitiveDetector(iname), m_isCalo(calo)
24 {
25  // for CMS hits
26  m_namesOfSD.push_back(iname);
27 
28  // Geant4 hit collection
29  collectionName.insert(iname);
30 
31  // register sensitive detector
32  G4SDManager * SDman = G4SDManager::GetSDMpointer();
33  SDman->AddNewDetector(this);
34 
35  const std::vector<std::string>& lvNames = clg.logicalNames(iname);
36  std::stringstream ss;
37  for (auto & lvname : lvNames) {
38  this->AssignSD(lvname);
39  ss << " " << lvname;
40  }
41  edm::LogVerbatim("SensitiveDetector")
42  << " <" << iname <<"> : Assigns SD to LVs " << ss.str();
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  const G4TouchableHistory * theTouchable=static_cast<const 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  const G4ThreeVector localCoordinates = preStepPoint->GetTouchable()->GetHistory()
80  ->GetTopTransform().TransformPoint(globalCoordinates);
81  return ConvertToLocal3DPoint(localCoordinates);
82 }
83 
85 {
86  const G4StepPoint * preStepPoint = step->GetPreStepPoint();
87  G4ThreeVector localCoordinates = preStepPoint->GetTouchable()->GetHistory()
88  ->GetTopTransform().TransformPoint(preStepPoint->GetPosition());
89  return ConvertToLocal3DPoint(localCoordinates);
90 }
91 
93 {
94  const G4ThreeVector& globalCoordinates = step->GetPostStepPoint()->GetPosition();
95  G4ThreeVector localCoordinates = step->GetPreStepPoint()->GetTouchable()->GetHistory()
96  ->GetTopTransform().TransformPoint(globalCoordinates);
97  return ConvertToLocal3DPoint(localCoordinates);
98 }
99 
101 {
102  TrackInformation* info = (TrackInformation*)(aTrack->GetUserInformation());
103  if (!info) {
104  edm::LogWarning("SensitiveDetector")
105  << " no TrackInformation available for trackID= "
106  << aTrack->GetTrackID();
107  throw SimG4Exception("SimG4CoreSensitiveDetector: cannot handle hits for "
108  + GetName());
109  }
110  return info;
111 }
112 
113 void SensitiveDetector::setNames(const std::vector<std::string>& hnames)
114 {
115  m_namesOfSD.clear();
116  m_namesOfSD = hnames;
117 }
118 
119 void SensitiveDetector::NaNTrap(const G4Step* aStep) const
120 {
121  const G4Track* CurrentTrk = aStep->GetTrack();
122  const G4ThreeVector& CurrentPos = CurrentTrk->GetPosition();
123  double xyz = CurrentPos.x() + CurrentPos.y() + CurrentPos.z();
124  const G4ThreeVector& CurrentMom = CurrentTrk->GetMomentum();
125  xyz += CurrentMom.x() + CurrentMom.y() + CurrentMom.z();
126 
127  if( edm::isNotFinite(xyz) ) {
128 
129  const G4VPhysicalVolume* pCurrentVol = aStep->GetPreStepPoint()->GetPhysicalVolume();
130  const G4String& NameOfVol = pCurrentVol->GetName();
131  throw SimG4Exception("SimG4CoreSensitiveDetector: Corrupted Event - NaN detected in volume "
132  + NameOfVol);
133  }
134 }
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 > & logicalNames(const std::string &readoutName) const
Compact representation of the geometrical detector hierarchy.
Definition: DDCompactView.h:80
void AssignSD(const std::string &vname)
std::vector< std::string > m_namesOfSD
std::string const collectionName[nCollections]
Definition: Collections.h:47
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 DDCompactView &cpv, 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