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 SensitiveDetectorCatalog &, 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 19 of file SensitiveDetector.h.

Member Enumeration Documentation

◆ coordinates

Enumerator
WorldCoordinates 
LocalCoordinates 

Definition at line 38 of file SensitiveDetector.h.

Constructor & Destructor Documentation

◆ SensitiveDetector()

SensitiveDetector::SensitiveDetector ( const std::string &  iname,
const SensitiveDetectorCatalog clg,
bool  calo 
)
explicit

Definition at line 17 of file SensitiveDetector.cc.

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

18  : G4VSensitiveDetector(iname), m_isCalo(calo) {
19  // for CMS hits
20  m_namesOfSD.push_back(iname);
21 
22  // Geant4 hit collection
23  collectionName.insert(iname);
24 
25  // register sensitive detector
26  G4SDManager* SDman = G4SDManager::GetSDMpointer();
27  SDman->AddNewDetector(this);
28 
29  const std::vector<std::string_view>& lvNames = clg.logicalNames(iname);
30  std::stringstream ss;
31  for (auto& lvname : lvNames) {
32  this->AssignSD({lvname.data(), lvname.size()});
33  ss << " " << lvname;
34  }
35  edm::LogVerbatim("SensitiveDetector") << " <" << iname << "> : Assigns SD to LVs " << ss.str();
36 }
Log< level::Info, true > LogVerbatim
const std::vector< std::string_view > logicalNames(const std::string &readoutName) const
void AssignSD(const std::string &vname)
std::vector< std::string > m_namesOfSD
Definition: Common.h:9

◆ ~SensitiveDetector()

SensitiveDetector::~SensitiveDetector ( )
override

Definition at line 38 of file SensitiveDetector.cc.

38 {}

Member Function Documentation

◆ AssignSD()

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

Definition at line 44 of file SensitiveDetector.cc.

Referenced by SensitiveDetector().

44  {
45  G4LogicalVolumeStore* theStore = G4LogicalVolumeStore::GetInstance();
46  for (auto& lv : *theStore) {
47  if (vname == lv->GetName()) {
48  lv->SetSensitiveDetector(this);
49  }
50  }
51 }

◆ clearHits()

virtual void SensitiveDetector::clearHits ( )
pure virtual

◆ cmsTrackInformation()

TrackInformation * SensitiveDetector::cmsTrackInformation ( const G4Track *  aTrack)
protected

Definition at line 91 of file SensitiveDetector.cc.

References info().

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

91  {
92  TrackInformation* info = (TrackInformation*)(aTrack->GetUserInformation());
93  if (nullptr == info) {
94  edm::LogWarning("SensitiveDetector") << " no TrackInformation available for trackID= " << aTrack->GetTrackID()
95  << " inside SD " << GetName();
96  G4Exception(
97  "SensitiveDetector::cmsTrackInformation()", "sd01", FatalException, "cannot handle hits without trackinfo");
98  }
99  return info;
100 }
static const TGPicture * info(bool iBackgroundIsBlack)
Log< level::Warning, false > LogWarning

◆ ConvertToLocal3DPoint()

Local3DPoint SensitiveDetector::ConvertToLocal3DPoint ( const G4ThreeVector &  point) const
inlineprotected

Definition at line 46 of file SensitiveDetector.h.

References point.

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

46  {
47  return Local3DPoint(point.x(), point.y(), point.z());
48  }
*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

◆ EndOfEvent()

void SensitiveDetector::EndOfEvent ( G4HCofThisEvent *  eventHC)
override

Definition at line 42 of file SensitiveDetector.cc.

42 {}

◆ FinalStepPosition()

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

Definition at line 65 of file SensitiveDetector.cc.

References hippyaddtobaddatafiles::cd(), ConvertToLocal3DPoint(), and WorldCoordinates.

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

65  {
66  const G4StepPoint* postStepPoint = step->GetPostStepPoint();
67  const G4ThreeVector& globalCoordinates = postStepPoint->GetPosition();
68  if (cd == WorldCoordinates) {
69  return ConvertToLocal3DPoint(globalCoordinates);
70  }
71  const G4StepPoint* preStepPoint = step->GetPreStepPoint();
72  const G4ThreeVector localCoordinates =
73  preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(globalCoordinates);
74  return ConvertToLocal3DPoint(localCoordinates);
75 }
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
step
Definition: StallMonitor.cc:98

◆ getNames()

const std::vector<std::string>& SensitiveDetector::getNames ( ) const
inline

Definition at line 32 of file SensitiveDetector.h.

References m_namesOfSD.

32 { return m_namesOfSD; }
std::vector< std::string > m_namesOfSD

◆ Initialize()

void SensitiveDetector::Initialize ( G4HCofThisEvent *  eventHC)
override

Definition at line 40 of file SensitiveDetector.cc.

40 {}

◆ InitialStepPosition()

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

Definition at line 53 of file SensitiveDetector.cc.

References hippyaddtobaddatafiles::cd(), ConvertToLocal3DPoint(), and WorldCoordinates.

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

53  {
54  const G4StepPoint* preStepPoint = step->GetPreStepPoint();
55  const G4ThreeVector& globalCoordinates = preStepPoint->GetPosition();
56  if (cd == WorldCoordinates) {
57  return ConvertToLocal3DPoint(globalCoordinates);
58  }
59  const G4TouchableHistory* theTouchable = static_cast<const G4TouchableHistory*>(preStepPoint->GetTouchable());
60  const G4ThreeVector localCoordinates =
61  theTouchable->GetHistory()->GetTopTransform().TransformPoint(globalCoordinates);
62  return ConvertToLocal3DPoint(localCoordinates);
63 }
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
step
Definition: StallMonitor.cc:98

◆ isCaloSD()

bool SensitiveDetector::isCaloSD ( ) const
inline

Definition at line 34 of file SensitiveDetector.h.

References m_isCalo.

34 { return m_isCalo; }

◆ LocalPostStepPosition()

Local3DPoint SensitiveDetector::LocalPostStepPosition ( const G4Step *  step) const
protected

Definition at line 84 of file SensitiveDetector.cc.

References ConvertToLocal3DPoint().

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

84  {
85  const G4ThreeVector& globalCoordinates = step->GetPostStepPoint()->GetPosition();
86  G4ThreeVector localCoordinates =
87  step->GetPreStepPoint()->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(globalCoordinates);
88  return ConvertToLocal3DPoint(localCoordinates);
89 }
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
step
Definition: StallMonitor.cc:98

◆ LocalPreStepPosition()

Local3DPoint SensitiveDetector::LocalPreStepPosition ( const G4Step *  step) const
protected

Definition at line 77 of file SensitiveDetector.cc.

References ConvertToLocal3DPoint().

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

77  {
78  const G4StepPoint* preStepPoint = step->GetPreStepPoint();
79  G4ThreeVector localCoordinates =
80  preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(preStepPoint->GetPosition());
81  return ConvertToLocal3DPoint(localCoordinates);
82 }
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
step
Definition: StallMonitor.cc:98

◆ NaNTrap()

void SensitiveDetector::NaNTrap ( const G4Step *  step) const
protected

Definition at line 107 of file SensitiveDetector.cc.

References edm::isNotFinite().

Referenced by CaloSD::ProcessHits().

107  {
108  G4Track* currentTrk = aStep->GetTrack();
109  double ekin = currentTrk->GetKineticEnergy();
110  if (ekin < 0.0) {
111  const G4VPhysicalVolume* pCurrentVol = aStep->GetPreStepPoint()->GetPhysicalVolume();
112  edm::LogWarning("SensitiveDetector") << "Negative kinetic energy Ekin(MeV)=" << ekin / CLHEP::MeV << " of "
113  << currentTrk->GetDefinition()->GetParticleName()
114  << " trackID= " << currentTrk->GetTrackID() << " inside "
115  << pCurrentVol->GetName();
116  currentTrk->SetKineticEnergy(0.0);
117  }
118  const G4ThreeVector& currentPos = currentTrk->GetPosition();
119  double xyz = currentPos.x() + currentPos.y() + currentPos.z();
120  const G4ThreeVector& currentMom = currentTrk->GetMomentum();
121  xyz += currentMom.x() + currentMom.y() + currentMom.z();
122 
123  if (edm::isNotFinite(xyz)) {
124  const G4VPhysicalVolume* pCurrentVol = aStep->GetPreStepPoint()->GetPhysicalVolume();
125  edm::LogWarning("SensitiveDetector") << "NaN detected for trackID= " << currentTrk->GetTrackID() << " inside "
126  << pCurrentVol->GetName();
127  G4Exception("SensitiveDetector::NaNTrap()", "sd01", FatalException, "corrupted event or step");
128  }
129 }
constexpr bool isNotFinite(T x)
Definition: isFinite.h:9
Log< level::Warning, false > LogWarning

◆ ProcessHits()

G4bool SensitiveDetector::ProcessHits ( G4Step *  step,
G4TouchableHistory *  tHistory 
)
overridepure virtual

◆ setDetUnitId()

virtual uint32_t SensitiveDetector::setDetUnitId ( const G4Step *  step)
pure virtual

◆ setNames()

void SensitiveDetector::setNames ( const std::vector< std::string > &  hnames)
protected

Definition at line 102 of file SensitiveDetector.cc.

References m_namesOfSD.

Referenced by TkAccumulatingSensitiveDetector::TkAccumulatingSensitiveDetector().

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

Member Data Documentation

◆ m_isCalo

bool SensitiveDetector::m_isCalo
private

Definition at line 59 of file SensitiveDetector.h.

Referenced by isCaloSD().

◆ m_namesOfSD

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

Definition at line 58 of file SensitiveDetector.h.

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