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  : G4VSensitiveDetector(iname), m_isCalo(calo) {
20  // for CMS hits
21  m_namesOfSD.push_back(iname);
22 
23  // Geant4 hit collection
24  collectionName.insert(iname);
25 
26  // register sensitive detector
27  G4SDManager* SDman = G4SDManager::GetSDMpointer();
28  SDman->AddNewDetector(this);
29 
30  const std::vector<std::string_view>& lvNames = clg.logicalNames(iname);
31  std::stringstream ss;
32  for (auto& lvname : lvNames) {
33  this->AssignSD({lvname.data(), lvname.size()});
34  ss << " " << lvname;
35  }
36  edm::LogVerbatim("SensitiveDetector") << " <" << iname << "> : Assigns SD to LVs " << ss.str();
37 }
38 
40 
41 void SensitiveDetector::Initialize(G4HCofThisEvent* eventHC) {}
42 
43 void SensitiveDetector::EndOfEvent(G4HCofThisEvent* eventHC) {}
44 
46  G4LogicalVolumeStore* theStore = G4LogicalVolumeStore::GetInstance();
47  for (auto& lv : *theStore) {
48  if (vname == lv->GetName()) {
49  lv->SetSensitiveDetector(this);
50  }
51  }
52 }
53 
55  const G4StepPoint* preStepPoint = step->GetPreStepPoint();
56  const G4ThreeVector& globalCoordinates = preStepPoint->GetPosition();
57  if (cd == WorldCoordinates) {
58  return ConvertToLocal3DPoint(globalCoordinates);
59  }
60  const G4TouchableHistory* theTouchable = static_cast<const G4TouchableHistory*>(preStepPoint->GetTouchable());
61  const G4ThreeVector localCoordinates =
62  theTouchable->GetHistory()->GetTopTransform().TransformPoint(globalCoordinates);
63  return ConvertToLocal3DPoint(localCoordinates);
64 }
65 
67  const G4StepPoint* postStepPoint = step->GetPostStepPoint();
68  const G4ThreeVector& globalCoordinates = postStepPoint->GetPosition();
69  if (cd == WorldCoordinates) {
70  return ConvertToLocal3DPoint(globalCoordinates);
71  }
72  const G4StepPoint* preStepPoint = step->GetPreStepPoint();
73  const G4ThreeVector localCoordinates =
74  preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(globalCoordinates);
75  return ConvertToLocal3DPoint(localCoordinates);
76 }
77 
79  const G4StepPoint* preStepPoint = step->GetPreStepPoint();
80  G4ThreeVector localCoordinates =
81  preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(preStepPoint->GetPosition());
82  return ConvertToLocal3DPoint(localCoordinates);
83 }
84 
86  const G4ThreeVector& globalCoordinates = step->GetPostStepPoint()->GetPosition();
87  G4ThreeVector localCoordinates =
88  step->GetPreStepPoint()->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(globalCoordinates);
89  return ConvertToLocal3DPoint(localCoordinates);
90 }
91 
93  TrackInformation* info = (TrackInformation*)(aTrack->GetUserInformation());
94  if (!info) {
95  edm::LogWarning("SensitiveDetector") << " no TrackInformation available for trackID= " << aTrack->GetTrackID();
96  throw SimG4Exception("SimG4CoreSensitiveDetector: cannot handle hits for " + GetName());
97  }
98  return info;
99 }
100 
101 void SensitiveDetector::setNames(const std::vector<std::string>& hnames) {
102  m_namesOfSD.clear();
103  m_namesOfSD = hnames;
104 }
105 
106 void SensitiveDetector::NaNTrap(const G4Step* aStep) const {
107  const G4Track* CurrentTrk = aStep->GetTrack();
108  const G4ThreeVector& CurrentPos = CurrentTrk->GetPosition();
109  double xyz = CurrentPos.x() + CurrentPos.y() + CurrentPos.z();
110  const G4ThreeVector& CurrentMom = CurrentTrk->GetMomentum();
111  xyz += CurrentMom.x() + CurrentMom.y() + CurrentMom.z();
112 
113  if (edm::isNotFinite(xyz)) {
114  const G4VPhysicalVolume* pCurrentVol = aStep->GetPreStepPoint()->GetPhysicalVolume();
115  const G4String& NameOfVol = pCurrentVol->GetName();
116  throw SimG4Exception("SimG4CoreSensitiveDetector: Corrupted Event - NaN detected in volume " + NameOfVol);
117  }
118 }
SensitiveDetector::SensitiveDetector
SensitiveDetector(const std::string &iname, const SensitiveDetectorCatalog &, bool calo)
Definition: SensitiveDetector.cc:18
MessageLogger.h
SensitiveDetector::~SensitiveDetector
~SensitiveDetector() override
Definition: SensitiveDetector.cc:39
step
step
Definition: StallMonitor.cc:94
SensitiveDetector::NaNTrap
void NaNTrap(const G4Step *step) const
Definition: SensitiveDetector.cc:106
edm::isNotFinite
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
materialBudgetVolume_cfi.lvNames
lvNames
Definition: materialBudgetVolume_cfi.py:6
SensitiveDetector::AssignSD
void AssignSD(const std::string &vname)
Definition: SensitiveDetector.cc:45
SensitiveDetector::coordinates
coordinates
Definition: SensitiveDetector.h:38
SimG4Exception
Definition: SimG4Exception.h:13
info
static const TGPicture * info(bool iBackgroundIsBlack)
Definition: FWCollectionSummaryWidget.cc:153
edm::LogWarning
Log< level::Warning, false > LogWarning
Definition: MessageLogger.h:122
SensitiveDetector::m_namesOfSD
std::vector< std::string > m_namesOfSD
Definition: SensitiveDetector.h:58
contentValuesCheck.ss
ss
Definition: contentValuesCheck.py:33
SensitiveDetectorCatalog
Definition: SensitiveDetectorCatalog.h:10
Point3DBase< float, LocalTag >
bysipixelclustmulteventfilter_cfi.collectionName
collectionName
Definition: bysipixelclustmulteventfilter_cfi.py:5
SensitiveDetector::FinalStepPosition
Local3DPoint FinalStepPosition(const G4Step *step, coordinates) const
Definition: SensitiveDetector.cc:66
SensitiveDetector::LocalPreStepPosition
Local3DPoint LocalPreStepPosition(const G4Step *step) const
Definition: SensitiveDetector.cc:78
SensitiveDetector::InitialStepPosition
Local3DPoint InitialStepPosition(const G4Step *step, coordinates) const
Definition: SensitiveDetector.cc:54
SimG4Exception.h
SensitiveDetector::ConvertToLocal3DPoint
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
Definition: SensitiveDetector.h:46
TrackInformation
Definition: TrackInformation.h:12
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
SensitiveDetector::setNames
void setNames(const std::vector< std::string > &)
Definition: SensitiveDetector.cc:101
SensitiveDetector::cmsTrackInformation
TrackInformation * cmsTrackInformation(const G4Track *aTrack)
Definition: SensitiveDetector.cc:92
SensitiveDetector::Initialize
void Initialize(G4HCofThisEvent *eventHC) override
Definition: SensitiveDetector.cc:41
isFinite.h
SensitiveDetector.h
edm::LogVerbatim
Log< level::Info, true > LogVerbatim
Definition: MessageLogger.h:128
SensitiveDetector::LocalPostStepPosition
Local3DPoint LocalPostStepPosition(const G4Step *step) const
Definition: SensitiveDetector.cc:85
SensitiveDetector::WorldCoordinates
Definition: SensitiveDetector.h:38
SensitiveDetectorCatalog::logicalNames
const std::vector< std::string_view > logicalNames(const std::string &readoutName) const
Definition: SensitiveDetectorCatalog.cc:31
SensitiveDetector::EndOfEvent
void EndOfEvent(G4HCofThisEvent *eventHC) override
Definition: SensitiveDetector.cc:43
hippyaddtobaddatafiles.cd
def cd(newdir)
Definition: hippyaddtobaddatafiles.py:40
calo
Definition: Common.h:9