CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Types | Protected Member Functions | Private Member Functions | Private Attributes
SensitiveDetector Class Referenceabstract

#include <SensitiveDetector.h>

Inheritance diagram for SensitiveDetector:
SensitiveCaloDetector SensitiveTkDetector CaloSD CaloTrkProcessing FiberSD HFChamberSD HFWedgeSD Bcm1fSD BHMSD BscSD FastTimerSD FP420SD MuonSensitiveDetector PltSD TkAccumulatingSensitiveDetector TotemSD

Public Member Functions

virtual void clearHits ()=0
 
void EndOfEvent (G4HCofThisEvent *eventHC) override
 
const std::vector< std::string > & getNames () const
 
void Initialize (G4HCofThisEvent *eventHC) override
 
G4bool ProcessHits (G4Step *step, G4TouchableHistory *tHistory) override=0
 
 SensitiveDetector (const std::string &iname, const DDCompactView &cpv, const SensitiveDetectorCatalog &, edm::ParameterSet const &p)
 
virtual uint32_t setDetUnitId (const G4Step *step)=0
 
 ~SensitiveDetector () override
 

Protected Types

enum  coordinates { WorldCoordinates, LocalCoordinates }
 

Protected Member Functions

Local3DPoint ConvertToLocal3DPoint (const G4ThreeVector &point) const
 
Local3DPoint FinalStepPosition (const G4Step *step, coordinates) const
 
Local3DPoint InitialStepPosition (const G4Step *step, coordinates) const
 
Local3DPoint LocalPostStepPosition (const G4Step *step) const
 
Local3DPoint LocalPreStepPosition (const G4Step *step) const
 
void NaNTrap (const G4Step *step) const
 
void setNames (const std::vector< std::string > &)
 

Private Member Functions

void AssignSD (const std::string &vname)
 

Private Attributes

std::vector< std::string > namesOfSD
 

Detailed Description

Definition at line 20 of file SensitiveDetector.h.

Member Enumeration Documentation

Enumerator
WorldCoordinates 
LocalCoordinates 

Definition at line 40 of file SensitiveDetector.h.

Constructor & Destructor Documentation

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

Definition at line 16 of file SensitiveDetector.cc.

References AssignSD(), ecaldqm::collectionName, SensitiveDetectorCatalog::logicalNames(), and namesOfSD.

19  :
20  G4VSensitiveDetector(iname)
21 {
22  // for CMS hits
23  namesOfSD.push_back(iname);
24 
25  // Geant4 hit collection
26  collectionName.insert(iname);
27 
28  // register sensitive detector
29  G4SDManager * SDman = G4SDManager::GetSDMpointer();
30  SDman->AddNewDetector(this);
31 
32  const std::vector<std::string>& lvNames = clg.logicalNames(iname);
33  std::stringstream ss;
34  for (auto & lvname : lvNames) {
35  this->AssignSD(lvname);
36  ss << " " << lvname;
37  }
38  edm::LogInfo("SensitiveDetector") << " <" << iname <<"> : Assigns SD to LVs "
39  << ss.str();
40 
41 }
const std::vector< std::string > & logicalNames(const std::string &readoutName) const
void AssignSD(const std::string &vname)
std::string const collectionName[nCollections]
Definition: Collections.h:47
std::vector< std::string > namesOfSD
SensitiveDetector::~SensitiveDetector ( )
override

Definition at line 43 of file SensitiveDetector.cc.

43 {}

Member Function Documentation

void SensitiveDetector::AssignSD ( const std::string &  vname)
private

Definition at line 49 of file SensitiveDetector.cc.

Referenced by ConvertToLocal3DPoint(), and SensitiveDetector().

50 {
51  G4LogicalVolumeStore * theStore = G4LogicalVolumeStore::GetInstance();
52  for (auto & lv : *theStore)
53  {
54  if (vname == lv->GetName()) {
55  lv->SetSensitiveDetector(this);
56  }
57  }
58 }
virtual void SensitiveDetector::clearHits ( )
pure virtual
Local3DPoint SensitiveDetector::ConvertToLocal3DPoint ( const G4ThreeVector &  point) const
inlineprotected

Definition at line 47 of file SensitiveDetector.h.

References AssignSD(), NaNTrap(), setNames(), and AlCaHLTBitMon_QueryRunRegistry::string.

Referenced by TkAccumulatingSensitiveDetector::createHit(), PltSD::createHit(), Bcm1fSD::createHit(), MuonSensitiveDetector::createHit(), FinalStepPosition(), MuonSensitiveDetector::FinalStepPositionVsParent(), InitialStepPosition(), MuonSensitiveDetector::InitialStepPositionVsParent(), LocalPostStepPosition(), LocalPreStepPosition(), and TkAccumulatingSensitiveDetector::updateHit().

48  { return Local3DPoint(point.x(),point.y(),point.z()); }
Point3DBase< float, LocalTag > Local3DPoint
Definition: LocalPoint.h:9
*vegas h *****************************************************used in the default bin number in original ***version of VEGAS is ***a higher bin number might help to derive a more precise ***grade subtle point
Definition: invegas.h:5
void SensitiveDetector::EndOfEvent ( G4HCofThisEvent *  eventHC)
override

Definition at line 47 of file SensitiveDetector.cc.

47 {}
Local3DPoint SensitiveDetector::FinalStepPosition ( const G4Step *  step,
coordinates  cd 
) const
protected

Definition at line 71 of file SensitiveDetector.cc.

References ConvertToLocal3DPoint(), and WorldCoordinates.

Referenced by Bcm1fSD::createHit(), PltSD::createHit(), MuonSensitiveDetector::createHit(), Bcm1fSD::updateHit(), PltSD::updateHit(), and MuonSensitiveDetector::updateHit().

72 {
73  const G4StepPoint * postStepPoint = step->GetPostStepPoint();
74  const G4ThreeVector& globalCoordinates = postStepPoint->GetPosition();
75  if (cd == WorldCoordinates) { return ConvertToLocal3DPoint(globalCoordinates); }
76  const G4StepPoint * preStepPoint = step->GetPreStepPoint();
77  G4TouchableHistory * theTouchable = (G4TouchableHistory *)(preStepPoint->GetTouchable());
78  const G4ThreeVector localCoordinates = theTouchable->GetHistory()
79  ->GetTopTransform().TransformPoint(globalCoordinates);
80  return ConvertToLocal3DPoint(localCoordinates);
81 }
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
step
const std::vector<std::string>& SensitiveDetector::getNames ( ) const
inline
void SensitiveDetector::Initialize ( G4HCofThisEvent *  eventHC)
override

Definition at line 45 of file SensitiveDetector.cc.

45 {}
Local3DPoint SensitiveDetector::InitialStepPosition ( const G4Step *  step,
coordinates  cd 
) const
protected

Definition at line 60 of file SensitiveDetector.cc.

References ConvertToLocal3DPoint(), and WorldCoordinates.

Referenced by Bcm1fSD::closeHit(), PltSD::closeHit(), Bcm1fSD::createHit(), PltSD::createHit(), MuonSensitiveDetector::createHit(), and MuonSensitiveDetector::ProcessHits().

61 {
62  const G4StepPoint * preStepPoint = step->GetPreStepPoint();
63  const G4ThreeVector& globalCoordinates = preStepPoint->GetPosition();
64  if (cd == WorldCoordinates) { return ConvertToLocal3DPoint(globalCoordinates); }
65  G4TouchableHistory * theTouchable=(G4TouchableHistory *)(preStepPoint->GetTouchable());
66  const G4ThreeVector localCoordinates = theTouchable->GetHistory()
67  ->GetTopTransform().TransformPoint(globalCoordinates);
68  return ConvertToLocal3DPoint(localCoordinates);
69 }
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
step
Local3DPoint SensitiveDetector::LocalPostStepPosition ( const G4Step *  step) const
protected

Definition at line 92 of file SensitiveDetector.cc.

References ConvertToLocal3DPoint().

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

93 {
94  const G4ThreeVector& globalCoordinates = step->GetPostStepPoint()->GetPosition();
95  G4TouchableHistory * theTouchable = (G4TouchableHistory *)(step->GetPreStepPoint()->GetTouchable());
96  G4ThreeVector localCoordinates = theTouchable->GetHistory()
97  ->GetTopTransform().TransformPoint(globalCoordinates);
98  return ConvertToLocal3DPoint(localCoordinates);
99 }
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
step
Local3DPoint SensitiveDetector::LocalPreStepPosition ( const G4Step *  step) const
protected

Definition at line 83 of file SensitiveDetector.cc.

References ConvertToLocal3DPoint().

Referenced by TkAccumulatingSensitiveDetector::closeHit(), and TkAccumulatingSensitiveDetector::createHit().

84 {
85  const G4StepPoint * preStepPoint = step->GetPreStepPoint();
86  G4TouchableHistory * theTouchable=(G4TouchableHistory *)(preStepPoint->GetTouchable());
87  G4ThreeVector localCoordinates = theTouchable->GetHistory()
88  ->GetTopTransform().TransformPoint(preStepPoint->GetPosition());
89  return ConvertToLocal3DPoint(localCoordinates);
90 }
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
step
void SensitiveDetector::NaNTrap ( const G4Step *  step) const
protected

Definition at line 107 of file SensitiveDetector.cc.

References edm::isNotFinite().

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

108 {
109  const G4Track* CurrentTrk = aStep->GetTrack();
110  const G4ThreeVector& CurrentPos = CurrentTrk->GetPosition();
111  double xyz = CurrentPos.x() + CurrentPos.y() + CurrentPos.z();
112  const G4ThreeVector& CurrentMom = CurrentTrk->GetMomentum();
113  xyz += CurrentMom.x() + CurrentMom.y() + CurrentMom.z();
114 
115  if( edm::isNotFinite(xyz) ) {
116 
117  const G4VPhysicalVolume* pCurrentVol = aStep->GetPreStepPoint()->GetPhysicalVolume();
118  const G4String& NameOfVol = pCurrentVol->GetName();
119  throw SimG4Exception( "SimG4CoreSensitiveDetector: Corrupted Event - NaN detected in volume "
120  + NameOfVol);
121  }
122 }
bool isNotFinite(T x)
Definition: isFinite.h:10
G4bool SensitiveDetector::ProcessHits ( G4Step *  step,
G4TouchableHistory *  tHistory 
)
overridepure virtual
virtual uint32_t SensitiveDetector::setDetUnitId ( const G4Step *  step)
pure virtual
void SensitiveDetector::setNames ( const std::vector< std::string > &  hnames)
protected

Definition at line 101 of file SensitiveDetector.cc.

References namesOfSD.

Referenced by ConvertToLocal3DPoint(), and TkAccumulatingSensitiveDetector::TkAccumulatingSensitiveDetector().

102 {
103  namesOfSD.clear();
104  namesOfSD = hnames;
105 }
std::vector< std::string > namesOfSD

Member Data Documentation

std::vector<std::string> SensitiveDetector::namesOfSD
private

Definition at line 57 of file SensitiveDetector.h.

Referenced by getNames(), SensitiveDetector(), and setNames().