CMS 3D CMS Logo

SensitiveDetector.cc
Go to the documentation of this file.
2 
5 
6 #include "G4SDManager.hh"
7 #include "G4Step.hh"
8 #include "G4Track.hh"
9 #include "G4StepPoint.hh"
10 #include "G4Transform3D.hh"
11 #include "G4LogicalVolumeStore.hh"
12 #include "G4TouchableHistory.hh"
13 #include "G4VUserTrackInformation.hh"
14 
15 #include <sstream>
16 
18  const SensitiveDetectorCatalog& clg,
19  bool calo,
20  const std::string& newcollname)
21  : G4VSensitiveDetector(iname), m_isCalo(calo) {
22  // for CMS hits
23  m_namesOfSD.push_back(iname);
24 
25  // Geant4 hit collection
26  collectionName.insert(iname);
27  if (!newcollname.empty())
28  collectionName.insert(newcollname);
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 << "\n";
39  }
40  if (newcollname.empty())
41  ss << " Collection " << iname;
42  else
43  ss << " Collections " << iname << " and " << newcollname;
44  edm::LogVerbatim("SensitiveDetector") << " <" << iname << "> : Assigns SD to " << lvNames.size() << " LVs "
45  << ss.str();
46 }
47 
49 
50 void SensitiveDetector::Initialize(G4HCofThisEvent* eventHC) {}
51 
52 void SensitiveDetector::EndOfEvent(G4HCofThisEvent* eventHC) {}
53 
55  G4LogicalVolumeStore* theStore = G4LogicalVolumeStore::GetInstance();
56  for (auto& lv : *theStore) {
57  if (vname == lv->GetName()) {
58  lv->SetSensitiveDetector(this);
59  }
60  }
61 }
62 
64  const G4StepPoint* preStepPoint = step->GetPreStepPoint();
65  const G4ThreeVector& globalCoordinates = preStepPoint->GetPosition();
66  if (cd == WorldCoordinates) {
67  return ConvertToLocal3DPoint(globalCoordinates);
68  }
69  const G4TouchableHistory* theTouchable = static_cast<const G4TouchableHistory*>(preStepPoint->GetTouchable());
70  const G4ThreeVector localCoordinates =
71  theTouchable->GetHistory()->GetTopTransform().TransformPoint(globalCoordinates);
72  return ConvertToLocal3DPoint(localCoordinates);
73 }
74 
76  const G4StepPoint* postStepPoint = step->GetPostStepPoint();
77  const G4ThreeVector& globalCoordinates = postStepPoint->GetPosition();
78  if (cd == WorldCoordinates) {
79  return ConvertToLocal3DPoint(globalCoordinates);
80  }
81  const G4StepPoint* preStepPoint = step->GetPreStepPoint();
82  const G4ThreeVector localCoordinates =
83  preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(globalCoordinates);
84  return ConvertToLocal3DPoint(localCoordinates);
85 }
86 
88  const G4StepPoint* preStepPoint = step->GetPreStepPoint();
89  G4ThreeVector localCoordinates =
90  preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(preStepPoint->GetPosition());
91  return ConvertToLocal3DPoint(localCoordinates);
92 }
93 
95  const G4ThreeVector& globalCoordinates = step->GetPostStepPoint()->GetPosition();
96  G4ThreeVector localCoordinates =
97  step->GetPreStepPoint()->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(globalCoordinates);
98  return ConvertToLocal3DPoint(localCoordinates);
99 }
100 
102  TrackInformation* info = (TrackInformation*)(aTrack->GetUserInformation());
103  if (nullptr == info) {
104  edm::LogWarning("SensitiveDetector") << " no TrackInformation available for trackID= " << aTrack->GetTrackID()
105  << " inside SD " << GetName();
106  G4Exception(
107  "SensitiveDetector::cmsTrackInformation()", "sd01", FatalException, "cannot handle hits without trackinfo");
108  }
109  return info;
110 }
111 
112 void SensitiveDetector::setNames(const std::vector<std::string>& hnames) {
113  m_namesOfSD.clear();
114  m_namesOfSD = hnames;
115 }
116 
117 void SensitiveDetector::NaNTrap(const G4Step* aStep) const {
118  G4Track* currentTrk = aStep->GetTrack();
119  double ekin = currentTrk->GetKineticEnergy();
120  if (ekin < 0.0) {
121  const G4VPhysicalVolume* pCurrentVol = aStep->GetPreStepPoint()->GetPhysicalVolume();
122  edm::LogWarning("SensitiveDetector") << "Negative kinetic energy Ekin(MeV)=" << ekin / CLHEP::MeV << " of "
123  << currentTrk->GetDefinition()->GetParticleName()
124  << " trackID= " << currentTrk->GetTrackID() << " inside "
125  << pCurrentVol->GetName();
126  currentTrk->SetKineticEnergy(0.0);
127  }
128  const G4ThreeVector& currentPos = currentTrk->GetPosition();
129  double xyz = currentPos.x() + currentPos.y() + currentPos.z();
130  const G4ThreeVector& currentMom = currentTrk->GetMomentum();
131  xyz += currentMom.x() + currentMom.y() + currentMom.z();
132 
133  if (edm::isNotFinite(xyz)) {
134  const G4VPhysicalVolume* pCurrentVol = aStep->GetPreStepPoint()->GetPhysicalVolume();
135  edm::LogWarning("SensitiveDetector") << "NaN detected for trackID= " << currentTrk->GetTrackID() << " inside "
136  << pCurrentVol->GetName();
137  G4Exception("SensitiveDetector::NaNTrap()", "sd01", FatalException, "corrupted event or step");
138  }
139 }
Log< level::Info, true > LogVerbatim
static const TGPicture * info(bool iBackgroundIsBlack)
~SensitiveDetector() override
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
SensitiveDetector(const std::string &iname, const SensitiveDetectorCatalog &, bool calo, const std::string &newcollname="")
const std::vector< std::string_view > logicalNames(const std::string &readoutName) const
void AssignSD(const std::string &vname)
Local3DPoint FinalStepPosition(const G4Step *step, coordinates) const
std::vector< std::string > m_namesOfSD
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
Local3DPoint InitialStepPosition(const G4Step *step, coordinates) const
TrackInformation * cmsTrackInformation(const G4Track *aTrack)
Local3DPoint LocalPreStepPosition(const G4Step *step) const
void EndOfEvent(G4HCofThisEvent *eventHC) override
void setNames(const std::vector< std::string > &)
Local3DPoint LocalPostStepPosition(const G4Step *step) const
step
Definition: StallMonitor.cc:83
Log< level::Warning, false > LogWarning
void Initialize(G4HCofThisEvent *eventHC) override
void NaNTrap(const G4Step *step) const
Definition: Common.h:9