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::LogVerbatim("SensitiveDetector")
41  << " <" << iname <<"> : Assigns SD to LVs " << ss.str();
42 }
43 
45 
46 void SensitiveDetector::Initialize(G4HCofThisEvent * eventHC) {}
47 
48 void SensitiveDetector::EndOfEvent(G4HCofThisEvent * eventHC) {}
49 
51 {
52  G4LogicalVolumeStore * theStore = G4LogicalVolumeStore::GetInstance();
53  for (auto & lv : *theStore)
54  {
55  if (vname == lv->GetName()) {
56  lv->SetSensitiveDetector(this);
57  }
58  }
59 }
60 
62 {
63  const G4StepPoint * preStepPoint = step->GetPreStepPoint();
64  const G4ThreeVector& globalCoordinates = preStepPoint->GetPosition();
65  if (cd == WorldCoordinates) { return ConvertToLocal3DPoint(globalCoordinates); }
66  const G4TouchableHistory * theTouchable=static_cast<const G4TouchableHistory *>(preStepPoint->GetTouchable());
67  const G4ThreeVector localCoordinates = theTouchable->GetHistory()
68  ->GetTopTransform().TransformPoint(globalCoordinates);
69  return ConvertToLocal3DPoint(localCoordinates);
70 }
71 
73 {
74  const G4StepPoint * postStepPoint = step->GetPostStepPoint();
75  const G4ThreeVector& globalCoordinates = postStepPoint->GetPosition();
76  if (cd == WorldCoordinates) { return ConvertToLocal3DPoint(globalCoordinates); }
77  const G4StepPoint * preStepPoint = step->GetPreStepPoint();
78  const G4ThreeVector localCoordinates = preStepPoint->GetTouchable()->GetHistory()
79  ->GetTopTransform().TransformPoint(globalCoordinates);
80  return ConvertToLocal3DPoint(localCoordinates);
81 }
82 
84 {
85  const G4StepPoint * preStepPoint = step->GetPreStepPoint();
86  G4ThreeVector localCoordinates = preStepPoint->GetTouchable()->GetHistory()
87  ->GetTopTransform().TransformPoint(preStepPoint->GetPosition());
88  return ConvertToLocal3DPoint(localCoordinates);
89 }
90 
92 {
93  const G4ThreeVector& globalCoordinates = step->GetPostStepPoint()->GetPosition();
94  G4ThreeVector localCoordinates = step->GetPreStepPoint()->GetTouchable()->GetHistory()
95  ->GetTopTransform().TransformPoint(globalCoordinates);
96  return ConvertToLocal3DPoint(localCoordinates);
97 }
98 
100 {
101  TrackInformation* info = (TrackInformation*)(aTrack->GetUserInformation());
102  if (!info) {
103  edm::LogWarning("SensitiveDetector")
104  << " no TrackInformation available for trackID= "
105  << aTrack->GetTrackID();
106  throw SimG4Exception("SimG4CoreSensitiveDetector: cannot handle hits for "
107  + GetName());
108  }
109  return info;
110 }
111 
112 void SensitiveDetector::setNames(const std::vector<std::string>& hnames)
113 {
114  namesOfSD.clear();
115  namesOfSD = hnames;
116 }
117 
118 void SensitiveDetector::NaNTrap(const G4Step* aStep) const
119 {
120  const G4Track* CurrentTrk = aStep->GetTrack();
121  const G4ThreeVector& CurrentPos = CurrentTrk->GetPosition();
122  double xyz = CurrentPos.x() + CurrentPos.y() + CurrentPos.z();
123  const G4ThreeVector& CurrentMom = CurrentTrk->GetMomentum();
124  xyz += CurrentMom.x() + CurrentMom.y() + CurrentMom.z();
125 
126  if( edm::isNotFinite(xyz) ) {
127 
128  const G4VPhysicalVolume* pCurrentVol = aStep->GetPreStepPoint()->GetPhysicalVolume();
129  const G4String& NameOfVol = pCurrentVol->GetName();
130  throw SimG4Exception("SimG4CoreSensitiveDetector: Corrupted Event - NaN detected in volume "
131  + NameOfVol);
132  }
133 }
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:80
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