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 FP420SD MuonSensitiveDetector PPSDiamondSD PPSPixelSD TimingSD TkAccumulatingSensitiveDetector TotemRPSD 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
 
bool isCaloSD () const
 
G4bool ProcessHits (G4Step *step, G4TouchableHistory *tHistory) override=0
 
 SensitiveDetector (const std::string &iname, const edm::EventSetup &es, const SensitiveDetectorCatalog &, edm::ParameterSet const &p, bool calo)
 
virtual uint32_t setDetUnitId (const G4Step *step)=0
 
 ~SensitiveDetector () override
 

Protected Types

enum  coordinates { WorldCoordinates, LocalCoordinates }
 

Protected Member Functions

TrackInformationcmsTrackInformation (const G4Track *aTrack)
 
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

bool m_isCalo
 
std::vector< std::string > m_namesOfSD
 

Detailed Description

Definition at line 25 of file SensitiveDetector.h.

Member Enumeration Documentation

Enumerator
WorldCoordinates 
LocalCoordinates 

Definition at line 48 of file SensitiveDetector.h.

Constructor & Destructor Documentation

SensitiveDetector::SensitiveDetector ( const std::string &  iname,
const edm::EventSetup es,
const SensitiveDetectorCatalog clg,
edm::ParameterSet const &  p,
bool  calo 
)
explicit

Definition at line 18 of file SensitiveDetector.cc.

References AssignSD(), bysipixelclustmulteventfilter_cfi::collectionName, SensitiveDetectorCatalog::logicalNames(), m_namesOfSD, and contentValuesCheck::ss.

23  : G4VSensitiveDetector(iname), m_isCalo(calo) {
24  // for CMS hits
25  m_namesOfSD.push_back(iname);
26 
27  // Geant4 hit collection
28  collectionName.insert(iname);
29 
30  // register sensitive detector
31  G4SDManager* SDman = G4SDManager::GetSDMpointer();
32  SDman->AddNewDetector(this);
33 
34  const std::vector<std::string_view>& lvNames = clg.logicalNames(iname);
35  std::stringstream ss;
36  for (auto& lvname : lvNames) {
37  this->AssignSD({lvname.data(), lvname.size()});
38  ss << " " << lvname;
39  }
40  edm::LogVerbatim("SensitiveDetector") << " <" << iname << "> : Assigns SD to LVs " << ss.str();
41 }
const std::vector< std::string_view > logicalNames(const std::string &readoutName) const
void AssignSD(const std::string &vname)
std::vector< std::string > m_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 SensitiveDetector().

49  {
50  G4LogicalVolumeStore* theStore = G4LogicalVolumeStore::GetInstance();
51  for (auto& lv : *theStore) {
52  if (vname == lv->GetName()) {
53  lv->SetSensitiveDetector(this);
54  }
55  }
56 }
virtual void SensitiveDetector::clearHits ( )
pure virtual
TrackInformation * SensitiveDetector::cmsTrackInformation ( const G4Track *  aTrack)
protected

Definition at line 96 of file SensitiveDetector.cc.

References info().

Referenced by TkAccumulatingSensitiveDetector::createHit(), MuonSensitiveDetector::createHit(), CaloSD::createNewHit(), CaloSD::getResponseWt(), TimingSD::getStepInfo(), CaloSD::getTrackID(), CaloSD::setTrackID(), TkAccumulatingSensitiveDetector::update(), and CaloSD::update().

96  {
97  TrackInformation* info = (TrackInformation*)(aTrack->GetUserInformation());
98  if (!info) {
99  edm::LogWarning("SensitiveDetector") << " no TrackInformation available for trackID= " << aTrack->GetTrackID();
100  throw SimG4Exception("SimG4CoreSensitiveDetector: cannot handle hits for " + GetName());
101  }
102  return info;
103 }
static const TGPicture * info(bool iBackgroundIsBlack)
Local3DPoint SensitiveDetector::ConvertToLocal3DPoint ( const G4ThreeVector &  point) const
inlineprotected

Definition at line 56 of file SensitiveDetector.h.

References AlCaHLTBitMon_QueryRunRegistry::string.

Referenced by TkAccumulatingSensitiveDetector::createHit(), MuonSensitiveDetector::createHit(), TimingSD::EndOfEvent(), FinalStepPosition(), MuonSensitiveDetector::FinalStepPositionVsParent(), InitialStepPosition(), MuonSensitiveDetector::InitialStepPositionVsParent(), LocalPostStepPosition(), LocalPreStepPosition(), TotemRPSD::stepInfo(), PPSDiamondSD::stepInfo(), and TkAccumulatingSensitiveDetector::updateHit().

56  {
57  return Local3DPoint(point.x(), point.y(), point.z());
58  }
*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
Point3DBase< float, LocalTag > Local3DPoint
Definition: LocalPoint.h:9
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 70 of file SensitiveDetector.cc.

References ConvertToLocal3DPoint(), and WorldCoordinates.

Referenced by MuonSensitiveDetector::createHit(), PPSPixelSD::stepInfo(), and MuonSensitiveDetector::updateHit().

70  {
71  const G4StepPoint* postStepPoint = step->GetPostStepPoint();
72  const G4ThreeVector& globalCoordinates = postStepPoint->GetPosition();
73  if (cd == WorldCoordinates) {
74  return ConvertToLocal3DPoint(globalCoordinates);
75  }
76  const G4StepPoint* preStepPoint = step->GetPreStepPoint();
77  const G4ThreeVector localCoordinates =
78  preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(globalCoordinates);
79  return ConvertToLocal3DPoint(localCoordinates);
80 }
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
step
Definition: StallMonitor.cc:94
const std::vector<std::string>& SensitiveDetector::getNames ( ) const
inline

Definition at line 42 of file SensitiveDetector.h.

Referenced by HcalTB06BeamSD::HcalTB06BeamSD().

42 { return m_namesOfSD; }
std::vector< std::string > m_namesOfSD
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 58 of file SensitiveDetector.cc.

References ConvertToLocal3DPoint(), and WorldCoordinates.

Referenced by MuonSensitiveDetector::createHit(), MuonSensitiveDetector::ProcessHits(), and PPSPixelSD::stepInfo().

58  {
59  const G4StepPoint* preStepPoint = step->GetPreStepPoint();
60  const G4ThreeVector& globalCoordinates = preStepPoint->GetPosition();
61  if (cd == WorldCoordinates) {
62  return ConvertToLocal3DPoint(globalCoordinates);
63  }
64  const G4TouchableHistory* theTouchable = static_cast<const G4TouchableHistory*>(preStepPoint->GetTouchable());
65  const G4ThreeVector localCoordinates =
66  theTouchable->GetHistory()->GetTopTransform().TransformPoint(globalCoordinates);
67  return ConvertToLocal3DPoint(localCoordinates);
68 }
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
step
Definition: StallMonitor.cc:94
bool SensitiveDetector::isCaloSD ( ) const
inline

Definition at line 44 of file SensitiveDetector.h.

44 { return m_isCalo; }
Local3DPoint SensitiveDetector::LocalPostStepPosition ( const G4Step *  step) const
protected

Definition at line 89 of file SensitiveDetector.cc.

References ConvertToLocal3DPoint().

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

89  {
90  const G4ThreeVector& globalCoordinates = step->GetPostStepPoint()->GetPosition();
91  G4ThreeVector localCoordinates =
92  step->GetPreStepPoint()->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(globalCoordinates);
93  return ConvertToLocal3DPoint(localCoordinates);
94 }
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
step
Definition: StallMonitor.cc:94
Local3DPoint SensitiveDetector::LocalPreStepPosition ( const G4Step *  step) const
protected

Definition at line 82 of file SensitiveDetector.cc.

References ConvertToLocal3DPoint().

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

82  {
83  const G4StepPoint* preStepPoint = step->GetPreStepPoint();
84  G4ThreeVector localCoordinates =
85  preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(preStepPoint->GetPosition());
86  return ConvertToLocal3DPoint(localCoordinates);
87 }
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
step
Definition: StallMonitor.cc:94
void SensitiveDetector::NaNTrap ( const G4Step *  step) const
protected

Definition at line 110 of file SensitiveDetector.cc.

References edm::isNotFinite().

Referenced by CaloSD::ProcessHits().

110  {
111  const G4Track* CurrentTrk = aStep->GetTrack();
112  const G4ThreeVector& CurrentPos = CurrentTrk->GetPosition();
113  double xyz = CurrentPos.x() + CurrentPos.y() + CurrentPos.z();
114  const G4ThreeVector& CurrentMom = CurrentTrk->GetMomentum();
115  xyz += CurrentMom.x() + CurrentMom.y() + CurrentMom.z();
116 
117  if (edm::isNotFinite(xyz)) {
118  const G4VPhysicalVolume* pCurrentVol = aStep->GetPreStepPoint()->GetPhysicalVolume();
119  const G4String& NameOfVol = pCurrentVol->GetName();
120  throw SimG4Exception("SimG4CoreSensitiveDetector: Corrupted Event - NaN detected in volume " + NameOfVol);
121  }
122 }
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
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 105 of file SensitiveDetector.cc.

References m_namesOfSD.

Referenced by TkAccumulatingSensitiveDetector::TkAccumulatingSensitiveDetector().

105  {
106  m_namesOfSD.clear();
107  m_namesOfSD = hnames;
108 }
std::vector< std::string > m_namesOfSD

Member Data Documentation

bool SensitiveDetector::m_isCalo
private

Definition at line 69 of file SensitiveDetector.h.

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

Definition at line 68 of file SensitiveDetector.h.

Referenced by SensitiveDetector(), and setNames().