CMS 3D CMS Logo

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