CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_7/src/SimG4Core/SensitiveDetector/src/SensitiveDetector.cc

Go to the documentation of this file.
00001 #include "SimG4Core/SensitiveDetector/interface/SensitiveDetector.h"
00002 
00003 #include "G4SDManager.hh"
00004 #include "G4Step.hh"
00005 #include "G4StepPoint.hh"
00006 #include "G4Transform3D.hh"
00007 #include "G4LogicalVolumeStore.hh"
00008 
00009 #include "SimG4Core/Notification/interface/SimG4Exception.h"
00010 #include "FWCore/Utilities/interface/isFinite.h"
00011 
00012 using std::string;
00013 
00014 SensitiveDetector::SensitiveDetector(string & iname, const DDCompactView & cpv,
00015                                      SensitiveDetectorCatalog & clg, 
00016                                      edm::ParameterSet const & p) :
00017   G4VSensitiveDetector(iname), name(iname) {}
00018 
00019 SensitiveDetector::~SensitiveDetector() {}
00020 
00021 void SensitiveDetector::Initialize(G4HCofThisEvent * eventHC) {}
00022 
00023 void SensitiveDetector::Register()
00024 {
00025   G4SDManager * SDman = G4SDManager::GetSDMpointer();
00026   SDman->AddNewDetector(this);
00027 }
00028 
00029 void SensitiveDetector::AssignSD(string & vname)
00030 {
00031     G4LogicalVolumeStore * theStore = G4LogicalVolumeStore::GetInstance();
00032     G4LogicalVolumeStore::const_iterator it;
00033     for (it = theStore->begin(); it != theStore->end(); it++)
00034     {
00035         G4LogicalVolume * v = *it;
00036         if (vname==v->GetName()) v->SetSensitiveDetector(this);
00037     }
00038 }
00039 
00040 void SensitiveDetector::EndOfEvent(G4HCofThisEvent * eventHC) {}
00041 
00042 #include "G4TouchableHistory.hh"
00043 
00044 Local3DPoint SensitiveDetector::InitialStepPosition(G4Step * s, coordinates c)
00045 {
00046     currentStep = s;
00047     G4StepPoint * preStepPoint = currentStep->GetPreStepPoint();
00048     G4ThreeVector globalCoordinates = preStepPoint->GetPosition();
00049     if (c == WorldCoordinates) return ConvertToLocal3DPoint(globalCoordinates);
00050     G4TouchableHistory * theTouchable=(G4TouchableHistory *)
00051                                       (preStepPoint->GetTouchable());
00052     G4ThreeVector localCoordinates = theTouchable->GetHistory()
00053                   ->GetTopTransform().TransformPoint(globalCoordinates);
00054     return ConvertToLocal3DPoint(localCoordinates); 
00055 }
00056 
00057 Local3DPoint SensitiveDetector::FinalStepPosition(G4Step * s, coordinates c)
00058 {
00059     currentStep = s;
00060     G4StepPoint * postStepPoint = currentStep->GetPostStepPoint();
00061     G4StepPoint * preStepPoint  = currentStep->GetPreStepPoint();
00062     G4ThreeVector globalCoordinates = postStepPoint->GetPosition();
00063     if (c == WorldCoordinates) return ConvertToLocal3DPoint(globalCoordinates);
00064     G4TouchableHistory * theTouchable = (G4TouchableHistory *)
00065                                         (preStepPoint->GetTouchable());
00066     G4ThreeVector localCoordinates = theTouchable->GetHistory()
00067                   ->GetTopTransform().TransformPoint(globalCoordinates);
00068     return ConvertToLocal3DPoint(localCoordinates); 
00069 }
00070 
00071 Local3DPoint SensitiveDetector::ConvertToLocal3DPoint(G4ThreeVector p)
00072 {
00073     return Local3DPoint(p.x(),p.y(),p.z());
00074 }
00075 
00076 void SensitiveDetector::NaNTrap( G4Step* aStep )
00077 {
00078 
00079     if ( aStep == NULL ) return ;
00080     
00081     G4Track* CurrentTrk = aStep->GetTrack() ;
00082     G4ThreeVector CurrentPos = CurrentTrk->GetPosition() ;
00083     G4ThreeVector CurrentMom = CurrentTrk->GetMomentum() ;
00084     G4VPhysicalVolume* pCurrentVol = CurrentTrk->GetVolume() ;
00085     G4String NameOfVol ;
00086     if ( pCurrentVol != NULL )
00087     {
00088        NameOfVol = pCurrentVol->GetName() ;
00089     }
00090     else
00091     {
00092        NameOfVol = "CorruptedVolumeInfo" ;
00093     }
00094     
00095     // for simplicity... maybe edm::isNotFinite() will work on the 3-vector directly...
00096     //
00097     double xyz[3] ;
00098     xyz[0] = CurrentPos.x() ;
00099     xyz[1] = CurrentPos.y() ;
00100     xyz[2] = CurrentPos.z() ;
00101     
00102     //
00103     // this is another trick to check on a NaN, maybe it's even CPU-faster...
00104     // but ler's stick to system function edm::isNotFinite(...) for now
00105     //
00106     // if ( !(xyz[0]==xyz[0]) || !(xyz[1]==xyz[1]) || !(xyz[2]==xyz[2]) )
00107     if( edm::isNotFinite(xyz[0]+xyz[1]+xyz[2]) != 0 )
00108     {
00109        // std::cout << " NaN detected in volume " << NameOfVol << std::endl ;
00110        throw SimG4Exception( "SimG4CoreSensitiveDetector: Corrupted Event - NaN detected (position)" ) ;
00111     }
00112 
00113     xyz[0] = CurrentMom.x() ;
00114     xyz[1] = CurrentMom.y() ;
00115     xyz[2] = CurrentMom.z() ;
00116     if ( !(xyz[0]==xyz[0]) || !(xyz[1]==xyz[1]) || !(xyz[2]==xyz[2]) ||
00117          edm::isNotFinite(xyz[0]) != 0 || edm::isNotFinite(xyz[1]) != 0 || edm::isNotFinite(xyz[2]) != 0 )
00118     {
00119        std::cout << " NaN detected in volume " << NameOfVol << std::endl ;
00120        throw SimG4Exception( "SimG4CoreSensitiveDetector: Corrupted Event - NaN detected (3-momentum)" ) ;
00121     }
00122 
00123    return ;
00124 
00125 }
00126 
00127 
00128