CMS 3D CMS Logo

SensitiveDetector.cc
Go to the documentation of this file.
2 
3 #include "G4SDManager.hh"
4 #include "G4Step.hh"
5 #include "G4StepPoint.hh"
6 #include "G4Transform3D.hh"
7 #include "G4LogicalVolumeStore.hh"
8 
11 
12 using std::string;
13 
15  const DDCompactView & cpv,
16  const SensitiveDetectorCatalog & clg,
17  edm::ParameterSet const & p) :
18  G4VSensitiveDetector(iname), name(iname) {}
19 
21 
22 void SensitiveDetector::Initialize(G4HCofThisEvent * eventHC) {}
23 
25 {
26  G4SDManager * SDman = G4SDManager::GetSDMpointer();
27  SDman->AddNewDetector(this);
28 }
29 
31 {
32  G4LogicalVolumeStore * theStore = G4LogicalVolumeStore::GetInstance();
33  G4LogicalVolumeStore::const_iterator it;
34  for (it = theStore->begin(); it != theStore->end(); it++)
35  {
36  G4LogicalVolume * v = *it;
37  if (vname==v->GetName()) { v->SetSensitiveDetector(this); }
38  }
39 }
40 
41 void SensitiveDetector::EndOfEvent(G4HCofThisEvent * eventHC) {}
42 
43 #include "G4TouchableHistory.hh"
44 
46 {
47  currentStep = s;
48  G4StepPoint * preStepPoint = currentStep->GetPreStepPoint();
49  G4ThreeVector globalCoordinates = preStepPoint->GetPosition();
50  if (c == WorldCoordinates) return ConvertToLocal3DPoint(globalCoordinates);
51  G4TouchableHistory * theTouchable=(G4TouchableHistory *)
52  (preStepPoint->GetTouchable());
53  G4ThreeVector localCoordinates = theTouchable->GetHistory()
54  ->GetTopTransform().TransformPoint(globalCoordinates);
55  return ConvertToLocal3DPoint(localCoordinates);
56 }
57 
59 {
60  currentStep = s;
61  G4StepPoint * postStepPoint = currentStep->GetPostStepPoint();
62  G4StepPoint * preStepPoint = currentStep->GetPreStepPoint();
63  G4ThreeVector globalCoordinates = postStepPoint->GetPosition();
64  if (c == WorldCoordinates) return ConvertToLocal3DPoint(globalCoordinates);
65  G4TouchableHistory * theTouchable = (G4TouchableHistory *)
66  (preStepPoint->GetTouchable());
67  G4ThreeVector localCoordinates = theTouchable->GetHistory()
68  ->GetTopTransform().TransformPoint(globalCoordinates);
69  return ConvertToLocal3DPoint(localCoordinates);
70 }
71 
73 {
74  return Local3DPoint(p.x(),p.y(),p.z());
75 }
76 
77 void SensitiveDetector::NaNTrap( G4Step* aStep )
78 {
79 
80  if ( aStep == nullptr ) return ;
81 
82  G4Track* CurrentTrk = aStep->GetTrack() ;
83  G4ThreeVector CurrentPos = CurrentTrk->GetPosition() ;
84  G4ThreeVector CurrentMom = CurrentTrk->GetMomentum() ;
85  G4VPhysicalVolume* pCurrentVol = CurrentTrk->GetVolume() ;
86  G4String NameOfVol ;
87  if ( pCurrentVol != nullptr )
88  {
89  NameOfVol = pCurrentVol->GetName() ;
90  }
91  else
92  {
93  NameOfVol = "CorruptedVolumeInfo" ;
94  }
95 
96  // for simplicity... maybe edm::isNotFinite() will work on the
97  // 3-vector directly...
98 
99  double xyz[3] ;
100  xyz[0] = CurrentPos.x() ;
101  xyz[1] = CurrentPos.y() ;
102  xyz[2] = CurrentPos.z() ;
103 
104  //
105  // this is another trick to check on a NaN, maybe it's even CPU-faster...
106  // but ler's stick to system function edm::isNotFinite(...) for now
107  //
108  // if ( !(xyz[0]==xyz[0]) || !(xyz[1]==xyz[1]) || !(xyz[2]==xyz[2]) )
109  if( edm::isNotFinite(xyz[0]+xyz[1]+xyz[2]) != 0 )
110  {
111  throw SimG4Exception( "SimG4CoreSensitiveDetector: Corrupted Event - NaN detected (position) in volume " + NameOfVol);
112  }
113 
114  xyz[0] = CurrentMom.x() ;
115  xyz[1] = CurrentMom.y() ;
116  xyz[2] = CurrentMom.z() ;
117  if ( !(xyz[0]==xyz[0]) || !(xyz[1]==xyz[1]) || !(xyz[2]==xyz[2]) ||
118  edm::isNotFinite(xyz[0]) != 0 || edm::isNotFinite(xyz[1]) != 0 ||
119  edm::isNotFinite(xyz[2]) != 0 )
120  {
121  throw SimG4Exception( "SimG4CoreSensitiveDetector: Corrupted Event - NaN detected (3-momentum) in volume " + NameOfVol);
122  }
123 
124  return;
125 
126 }
virtual ~SensitiveDetector()
type of data representation of DDCompactView
Definition: DDCompactView.h:90
SensitiveDetector(std::string &iname, const DDCompactView &cpv, const SensitiveDetectorCatalog &, edm::ParameterSet const &p)
virtual void AssignSD(const std::string &vname)
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point)
bool isNotFinite(T x)
Definition: isFinite.h:10
void NaNTrap(G4Step *step)
Point3DBase< float, LocalTag > Local3DPoint
Definition: LocalPoint.h:9
return(e1-e2)*(e1-e2)+dp *dp
virtual void EndOfEvent(G4HCofThisEvent *eventHC)
Local3DPoint FinalStepPosition(G4Step *s, coordinates)
Local3DPoint InitialStepPosition(G4Step *s, coordinates)
virtual void Initialize(G4HCofThisEvent *eventHC)