CMS 3D CMS Logo

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