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 "G4StepPoint.hh"
10 #include "G4Transform3D.hh"
11 #include "G4LogicalVolumeStore.hh"
12 #include "G4TouchableHistory.hh"
13 
14 #include <sstream>
15 
17  const DDCompactView & cpv,
18  const SensitiveDetectorCatalog & clg,
19  edm::ParameterSet const & p) :
20  G4VSensitiveDetector(iname)
21 {
22  // for CMS hits
23  namesOfSD.push_back(iname);
24 
25  // Geant4 hit collection
26  collectionName.insert(iname);
27 
28  // register sensitive detector
29  G4SDManager * SDman = G4SDManager::GetSDMpointer();
30  SDman->AddNewDetector(this);
31 
32  const std::vector<std::string>& lvNames = clg.logicalNames(iname);
33  std::stringstream ss;
34  for (auto & lvname : lvNames) {
35  this->AssignSD(lvname);
36  ss << " " << lvname;
37  }
38  edm::LogInfo("SensitiveDetector") << " <" << iname <<"> : Assigns SD to LVs "
39  << ss.str();
40 
41 }
42 
44 
45 void SensitiveDetector::Initialize(G4HCofThisEvent * eventHC) {}
46 
47 void SensitiveDetector::EndOfEvent(G4HCofThisEvent * eventHC) {}
48 
50 {
51  G4LogicalVolumeStore * theStore = G4LogicalVolumeStore::GetInstance();
52  for (auto & lv : *theStore)
53  {
54  if (vname == lv->GetName()) {
55  lv->SetSensitiveDetector(this);
56  }
57  }
58 }
59 
61 {
62  const G4StepPoint * preStepPoint = step->GetPreStepPoint();
63  const G4ThreeVector& globalCoordinates = preStepPoint->GetPosition();
64  if (cd == WorldCoordinates) { return ConvertToLocal3DPoint(globalCoordinates); }
65  G4TouchableHistory * theTouchable=(G4TouchableHistory *)(preStepPoint->GetTouchable());
66  const G4ThreeVector localCoordinates = theTouchable->GetHistory()
67  ->GetTopTransform().TransformPoint(globalCoordinates);
68  return ConvertToLocal3DPoint(localCoordinates);
69 }
70 
72 {
73  const G4StepPoint * postStepPoint = step->GetPostStepPoint();
74  const G4ThreeVector& globalCoordinates = postStepPoint->GetPosition();
75  if (cd == WorldCoordinates) { return ConvertToLocal3DPoint(globalCoordinates); }
76  const G4StepPoint * preStepPoint = step->GetPreStepPoint();
77  G4TouchableHistory * theTouchable = (G4TouchableHistory *)(preStepPoint->GetTouchable());
78  const G4ThreeVector localCoordinates = theTouchable->GetHistory()
79  ->GetTopTransform().TransformPoint(globalCoordinates);
80  return ConvertToLocal3DPoint(localCoordinates);
81 }
82 
84 {
85  const G4StepPoint * preStepPoint = step->GetPreStepPoint();
86  G4TouchableHistory * theTouchable=(G4TouchableHistory *)(preStepPoint->GetTouchable());
87  G4ThreeVector localCoordinates = theTouchable->GetHistory()
88  ->GetTopTransform().TransformPoint(preStepPoint->GetPosition());
89  return ConvertToLocal3DPoint(localCoordinates);
90 }
91 
93 {
94  const G4ThreeVector& globalCoordinates = step->GetPostStepPoint()->GetPosition();
95  G4TouchableHistory * theTouchable = (G4TouchableHistory *)(step->GetPreStepPoint()->GetTouchable());
96  G4ThreeVector localCoordinates = theTouchable->GetHistory()
97  ->GetTopTransform().TransformPoint(globalCoordinates);
98  return ConvertToLocal3DPoint(localCoordinates);
99 }
100 
101 void SensitiveDetector::setNames(const std::vector<std::string>& hnames)
102 {
103  namesOfSD.clear();
104  namesOfSD = hnames;
105 }
106 
107 void SensitiveDetector::NaNTrap(const G4Step* aStep) const
108 {
109  const G4Track* CurrentTrk = aStep->GetTrack();
110  const G4ThreeVector& CurrentPos = CurrentTrk->GetPosition();
111  double xyz = CurrentPos.x() + CurrentPos.y() + CurrentPos.z();
112  const G4ThreeVector& CurrentMom = CurrentTrk->GetMomentum();
113  xyz += CurrentMom.x() + CurrentMom.y() + CurrentMom.z();
114 
115  if( edm::isNotFinite(xyz) ) {
116 
117  const G4VPhysicalVolume* pCurrentVol = aStep->GetPreStepPoint()->GetPhysicalVolume();
118  const G4String& NameOfVol = pCurrentVol->GetName();
119  throw SimG4Exception( "SimG4CoreSensitiveDetector: Corrupted Event - NaN detected in volume "
120  + NameOfVol);
121  }
122 }
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
type of data representation of DDCompactView
Definition: DDCompactView.h:90
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
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