CMS 3D CMS Logo

Public Types | Public Member Functions | Private Attributes

SensitiveDetector Class Reference

#include <SensitiveDetector.h>

Inheritance diagram for SensitiveDetector:
SensitiveCaloDetector SensitiveTkDetector CaloSD CaloTrkProcessing FiberSD HFChamberSD HFWedgeSD BscSD FP420SD MuonSensitiveDetector PLTSensitiveDetector TkAccumulatingSensitiveDetector TotemSD

List of all members.

Public Types

enum  coordinates { WorldCoordinates, LocalCoordinates }

Public Member Functions

virtual void AssignSD (std::string &vname)
virtual void clearHits ()=0
Local3DPoint ConvertToLocal3DPoint (G4ThreeVector point)
virtual void EndOfEvent (G4HCofThisEvent *eventHC)
Local3DPoint FinalStepPosition (G4Step *s, coordinates)
virtual std::vector< std::string > getNames ()
virtual void Initialize (G4HCofThisEvent *eventHC)
Local3DPoint InitialStepPosition (G4Step *s, coordinates)
std::string nameOfSD ()
void NaNTrap (G4Step *step)
virtual G4bool ProcessHits (G4Step *step, G4TouchableHistory *tHistory)=0
void Register ()
 SensitiveDetector (std::string &iname, const DDCompactView &cpv, SensitiveDetectorCatalog &, edm::ParameterSet const &p)
virtual uint32_t setDetUnitId (G4Step *step)=0
virtual ~SensitiveDetector ()

Private Attributes

G4Step * currentStep
std::string name

Detailed Description

Definition at line 24 of file SensitiveDetector.h.


Member Enumeration Documentation

Enumerator:
WorldCoordinates 
LocalCoordinates 

Definition at line 38 of file SensitiveDetector.h.


Constructor & Destructor Documentation

SensitiveDetector::SensitiveDetector ( std::string &  iname,
const DDCompactView cpv,
SensitiveDetectorCatalog ,
edm::ParameterSet const &  p 
) [explicit]

Definition at line 14 of file SensitiveDetector.cc.

                                                                :
  G4VSensitiveDetector(iname), name(iname) {}
SensitiveDetector::~SensitiveDetector ( ) [virtual]

Definition at line 19 of file SensitiveDetector.cc.

{}

Member Function Documentation

void SensitiveDetector::AssignSD ( std::string &  vname) [virtual]

Definition at line 29 of file SensitiveDetector.cc.

References v.

Referenced by BscSD::BscSD(), CaloSD::CaloSD(), FiberSD::FiberSD(), FP420SD::FP420SD(), HFChamberSD::HFChamberSD(), HFWedgeSD::HFWedgeSD(), MuonSensitiveDetector::MuonSensitiveDetector(), PLTSensitiveDetector::PLTSensitiveDetector(), TkAccumulatingSensitiveDetector::TkAccumulatingSensitiveDetector(), and TotemSD::TotemSD().

{
    G4LogicalVolumeStore * theStore = G4LogicalVolumeStore::GetInstance();
    G4LogicalVolumeStore::const_iterator it;
    for (it = theStore->begin(); it != theStore->end(); it++)
    {
        G4LogicalVolume * v = *it;
        if (vname==v->GetName()) v->SetSensitiveDetector(this);
    }
}
virtual void SensitiveDetector::clearHits ( ) [pure virtual]
Local3DPoint SensitiveDetector::ConvertToLocal3DPoint ( G4ThreeVector  point)
void SensitiveDetector::EndOfEvent ( G4HCofThisEvent *  eventHC) [virtual]
Local3DPoint SensitiveDetector::FinalStepPosition ( G4Step *  s,
coordinates  c 
)

Definition at line 57 of file SensitiveDetector.cc.

References ConvertToLocal3DPoint(), currentStep, alignCSCRings::s, and WorldCoordinates.

Referenced by MuonSensitiveDetector::createHit(), PLTSensitiveDetector::createHit(), TkAccumulatingSensitiveDetector::createHit(), MuonSensitiveDetector::updateHit(), PLTSensitiveDetector::updateHit(), and TkAccumulatingSensitiveDetector::updateHit().

{
    currentStep = s;
    G4StepPoint * postStepPoint = currentStep->GetPostStepPoint();
    G4StepPoint * preStepPoint  = currentStep->GetPreStepPoint();
    G4ThreeVector globalCoordinates = postStepPoint->GetPosition();
    if (c == WorldCoordinates) return ConvertToLocal3DPoint(globalCoordinates);
    G4TouchableHistory * theTouchable = (G4TouchableHistory *)
                                        (preStepPoint->GetTouchable());
    G4ThreeVector localCoordinates = theTouchable->GetHistory()
                  ->GetTopTransform().TransformPoint(globalCoordinates);
    return ConvertToLocal3DPoint(localCoordinates); 
}
virtual std::vector<std::string> SensitiveDetector::getNames ( ) [inline, virtual]

Reimplemented in BscSD, FP420SD, MuonSensitiveDetector, and TkAccumulatingSensitiveDetector.

Definition at line 43 of file SensitiveDetector.h.

References nameOfSD(), and groupFilesInBlocks::temp.

Referenced by CaloTrkProcessing::CaloTrkProcessing(), HCalSD::HCalSD(), and HcalTB06BeamSD::HcalTB06BeamSD().

  {
    std::vector<std::string> temp;
    temp.push_back(nameOfSD());
    return temp;
  }
void SensitiveDetector::Initialize ( G4HCofThisEvent *  eventHC) [virtual]

Reimplemented in CaloSD, CaloTrkProcessing, BscSD, TotemSD, FP420SD, FiberSD, HFChamberSD, and HFWedgeSD.

Definition at line 21 of file SensitiveDetector.cc.

{}
Local3DPoint SensitiveDetector::InitialStepPosition ( G4Step *  s,
coordinates  c 
)

Definition at line 44 of file SensitiveDetector.cc.

References ConvertToLocal3DPoint(), currentStep, alignCSCRings::s, and WorldCoordinates.

Referenced by TkAccumulatingSensitiveDetector::closeHit(), PLTSensitiveDetector::closeHit(), MuonSensitiveDetector::createHit(), PLTSensitiveDetector::createHit(), TkAccumulatingSensitiveDetector::createHit(), and MuonSensitiveDetector::ProcessHits().

{
    currentStep = s;
    G4StepPoint * preStepPoint = currentStep->GetPreStepPoint();
    G4ThreeVector globalCoordinates = preStepPoint->GetPosition();
    if (c == WorldCoordinates) return ConvertToLocal3DPoint(globalCoordinates);
    G4TouchableHistory * theTouchable=(G4TouchableHistory *)
                                      (preStepPoint->GetTouchable());
    G4ThreeVector localCoordinates = theTouchable->GetHistory()
                  ->GetTopTransform().TransformPoint(globalCoordinates);
    return ConvertToLocal3DPoint(localCoordinates); 
}
std::string SensitiveDetector::nameOfSD ( ) [inline]

Definition at line 42 of file SensitiveDetector.h.

References name.

Referenced by getNames().

{ return name; }
void SensitiveDetector::NaNTrap ( G4Step *  step)

Definition at line 76 of file SensitiveDetector.cc.

References gather_cfg::cout, edm::isNotFinite(), and NULL.

Referenced by CaloSD::ProcessHits(), ZdcSD::ProcessHits(), and HCalSD::ProcessHits().

{

    if ( aStep == NULL ) return ;
    
    G4Track* CurrentTrk = aStep->GetTrack() ;
    G4ThreeVector CurrentPos = CurrentTrk->GetPosition() ;
    G4ThreeVector CurrentMom = CurrentTrk->GetMomentum() ;
    G4VPhysicalVolume* pCurrentVol = CurrentTrk->GetVolume() ;
    G4String NameOfVol ;
    if ( pCurrentVol != NULL )
    {
       NameOfVol = pCurrentVol->GetName() ;
    }
    else
    {
       NameOfVol = "CorruptedVolumeInfo" ;
    }
    
    // for simplicity... maybe edm::isNotFinite() will work on the 3-vector directly...
    //
    double xyz[3] ;
    xyz[0] = CurrentPos.x() ;
    xyz[1] = CurrentPos.y() ;
    xyz[2] = CurrentPos.z() ;
    
    //
    // this is another trick to check on a NaN, maybe it's even CPU-faster...
    // but ler's stick to system function edm::isNotFinite(...) for now
    //
    // if ( !(xyz[0]==xyz[0]) || !(xyz[1]==xyz[1]) || !(xyz[2]==xyz[2]) )
    if( edm::isNotFinite(xyz[0]+xyz[1]+xyz[2]) != 0 )
    {
       // std::cout << " NaN detected in volume " << NameOfVol << std::endl ;
       throw SimG4Exception( "SimG4CoreSensitiveDetector: Corrupted Event - NaN detected (position)" ) ;
    }

    xyz[0] = CurrentMom.x() ;
    xyz[1] = CurrentMom.y() ;
    xyz[2] = CurrentMom.z() ;
    if ( !(xyz[0]==xyz[0]) || !(xyz[1]==xyz[1]) || !(xyz[2]==xyz[2]) ||
         edm::isNotFinite(xyz[0]) != 0 || edm::isNotFinite(xyz[1]) != 0 || edm::isNotFinite(xyz[2]) != 0 )
    {
       std::cout << " NaN detected in volume " << NameOfVol << std::endl ;
       throw SimG4Exception( "SimG4CoreSensitiveDetector: Corrupted Event - NaN detected (3-momentum)" ) ;
    }

   return ;

}
virtual G4bool SensitiveDetector::ProcessHits ( G4Step *  step,
G4TouchableHistory *  tHistory 
) [pure virtual]
void SensitiveDetector::Register ( )
virtual uint32_t SensitiveDetector::setDetUnitId ( G4Step *  step) [pure virtual]

Member Data Documentation

G4Step* SensitiveDetector::currentStep [private]

Definition at line 54 of file SensitiveDetector.h.

Referenced by FinalStepPosition(), and InitialStepPosition().

std::string SensitiveDetector::name [private]