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, const std::string &newcollname="")
 
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 41 of file SensitiveDetector.h.

Constructor & Destructor Documentation

◆ SensitiveDetector()

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

Definition at line 17 of file SensitiveDetector.cc.

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

21  : G4VSensitiveDetector(iname), m_isCalo(calo) {
22  // for CMS hits
23  m_namesOfSD.push_back(iname);
24 
25  // Geant4 hit collection
26  collectionName.insert(iname);
27  if (!newcollname.empty())
28  collectionName.insert(newcollname);
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 << "\n";
39  }
40  if (newcollname.empty())
41  ss << " Collection " << iname;
42  else
43  ss << " Collections " << iname << " and " << newcollname;
44  edm::LogVerbatim("SensitiveDetector") << " <" << iname << "> : Assigns SD to " << lvNames.size() << " LVs "
45  << ss.str();
46 }
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 48 of file SensitiveDetector.cc.

48 {}

Member Function Documentation

◆ AssignSD()

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

Definition at line 54 of file SensitiveDetector.cc.

Referenced by SensitiveDetector().

54  {
55  G4LogicalVolumeStore* theStore = G4LogicalVolumeStore::GetInstance();
56  for (auto& lv : *theStore) {
57  if (vname == lv->GetName()) {
58  lv->SetSensitiveDetector(this);
59  }
60  }
61 }

◆ clearHits()

virtual void SensitiveDetector::clearHits ( )
pure virtual

◆ cmsTrackInformation()

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

Definition at line 101 of file SensitiveDetector.cc.

References info().

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

101  {
102  TrackInformation* info = (TrackInformation*)(aTrack->GetUserInformation());
103  if (nullptr == info) {
104  edm::LogWarning("SensitiveDetector") << " no TrackInformation available for trackID= " << aTrack->GetTrackID()
105  << " inside SD " << GetName();
106  G4Exception(
107  "SensitiveDetector::cmsTrackInformation()", "sd01", FatalException, "cannot handle hits without trackinfo");
108  }
109  return info;
110 }
static const TGPicture * info(bool iBackgroundIsBlack)
Log< level::Warning, false > LogWarning

◆ ConvertToLocal3DPoint()

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

Definition at line 49 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().

49  {
50  return Local3DPoint(point.x(), point.y(), point.z());
51  }
*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 52 of file SensitiveDetector.cc.

52 {}

◆ FinalStepPosition()

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

Definition at line 75 of file SensitiveDetector.cc.

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

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

75  {
76  const G4StepPoint* postStepPoint = step->GetPostStepPoint();
77  const G4ThreeVector& globalCoordinates = postStepPoint->GetPosition();
78  if (cd == WorldCoordinates) {
79  return ConvertToLocal3DPoint(globalCoordinates);
80  }
81  const G4StepPoint* preStepPoint = step->GetPreStepPoint();
82  const G4ThreeVector localCoordinates =
83  preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(globalCoordinates);
84  return ConvertToLocal3DPoint(localCoordinates);
85 }
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
step
Definition: StallMonitor.cc:83

◆ getNames()

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

Definition at line 35 of file SensitiveDetector.h.

References m_namesOfSD.

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

◆ Initialize()

void SensitiveDetector::Initialize ( G4HCofThisEvent *  eventHC)
override

Definition at line 50 of file SensitiveDetector.cc.

50 {}

◆ InitialStepPosition()

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

Definition at line 63 of file SensitiveDetector.cc.

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

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

63  {
64  const G4StepPoint* preStepPoint = step->GetPreStepPoint();
65  const G4ThreeVector& globalCoordinates = preStepPoint->GetPosition();
66  if (cd == WorldCoordinates) {
67  return ConvertToLocal3DPoint(globalCoordinates);
68  }
69  const G4TouchableHistory* theTouchable = static_cast<const G4TouchableHistory*>(preStepPoint->GetTouchable());
70  const G4ThreeVector localCoordinates =
71  theTouchable->GetHistory()->GetTopTransform().TransformPoint(globalCoordinates);
72  return ConvertToLocal3DPoint(localCoordinates);
73 }
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
step
Definition: StallMonitor.cc:83

◆ isCaloSD()

bool SensitiveDetector::isCaloSD ( ) const
inline

Definition at line 37 of file SensitiveDetector.h.

References m_isCalo.

37 { return m_isCalo; }

◆ LocalPostStepPosition()

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

Definition at line 94 of file SensitiveDetector.cc.

References ConvertToLocal3DPoint().

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

94  {
95  const G4ThreeVector& globalCoordinates = step->GetPostStepPoint()->GetPosition();
96  G4ThreeVector localCoordinates =
97  step->GetPreStepPoint()->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(globalCoordinates);
98  return ConvertToLocal3DPoint(localCoordinates);
99 }
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
step
Definition: StallMonitor.cc:83

◆ LocalPreStepPosition()

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

Definition at line 87 of file SensitiveDetector.cc.

References ConvertToLocal3DPoint().

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

87  {
88  const G4StepPoint* preStepPoint = step->GetPreStepPoint();
89  G4ThreeVector localCoordinates =
90  preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(preStepPoint->GetPosition());
91  return ConvertToLocal3DPoint(localCoordinates);
92 }
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
step
Definition: StallMonitor.cc:83

◆ NaNTrap()

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

Definition at line 117 of file SensitiveDetector.cc.

References edm::isNotFinite().

Referenced by CaloSD::ProcessHits().

117  {
118  G4Track* currentTrk = aStep->GetTrack();
119  double ekin = currentTrk->GetKineticEnergy();
120  if (ekin < 0.0) {
121  const G4VPhysicalVolume* pCurrentVol = aStep->GetPreStepPoint()->GetPhysicalVolume();
122  edm::LogWarning("SensitiveDetector") << "Negative kinetic energy Ekin(MeV)=" << ekin / CLHEP::MeV << " of "
123  << currentTrk->GetDefinition()->GetParticleName()
124  << " trackID= " << currentTrk->GetTrackID() << " inside "
125  << pCurrentVol->GetName();
126  currentTrk->SetKineticEnergy(0.0);
127  }
128  const G4ThreeVector& currentPos = currentTrk->GetPosition();
129  double xyz = currentPos.x() + currentPos.y() + currentPos.z();
130  const G4ThreeVector& currentMom = currentTrk->GetMomentum();
131  xyz += currentMom.x() + currentMom.y() + currentMom.z();
132 
133  if (edm::isNotFinite(xyz)) {
134  const G4VPhysicalVolume* pCurrentVol = aStep->GetPreStepPoint()->GetPhysicalVolume();
135  edm::LogWarning("SensitiveDetector") << "NaN detected for trackID= " << currentTrk->GetTrackID() << " inside "
136  << pCurrentVol->GetName();
137  G4Exception("SensitiveDetector::NaNTrap()", "sd01", FatalException, "corrupted event or step");
138  }
139 }
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 112 of file SensitiveDetector.cc.

References m_namesOfSD.

Referenced by TkAccumulatingSensitiveDetector::TkAccumulatingSensitiveDetector().

112  {
113  m_namesOfSD.clear();
114  m_namesOfSD = hnames;
115 }
std::vector< std::string > m_namesOfSD

Member Data Documentation

◆ m_isCalo

bool SensitiveDetector::m_isCalo
private

Definition at line 62 of file SensitiveDetector.h.

Referenced by isCaloSD().

◆ m_namesOfSD

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

Definition at line 61 of file SensitiveDetector.h.

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