CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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

Enumerator
WorldCoordinates 
LocalCoordinates 

Definition at line 38 of file SensitiveDetector.h.

Constructor & Destructor Documentation

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

Definition at line 18 of file SensitiveDetector.cc.

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

19  : G4VSensitiveDetector(iname), m_isCalo(calo) {
20  // for CMS hits
21  m_namesOfSD.push_back(iname);
22 
23  // Geant4 hit collection
24  collectionName.insert(iname);
25 
26  // register sensitive detector
27  G4SDManager* SDman = G4SDManager::GetSDMpointer();
28  SDman->AddNewDetector(this);
29 
30  const std::vector<std::string_view>& lvNames = clg.logicalNames(iname);
31  std::stringstream ss;
32  for (auto& lvname : lvNames) {
33  this->AssignSD({lvname.data(), lvname.size()});
34  ss << " " << lvname;
35  }
36  edm::LogVerbatim("SensitiveDetector") << " <" << iname << "> : Assigns SD to LVs " << ss.str();
37 }
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
std::string const collectionName[nCollections]
Definition: Collections.h:47
SensitiveDetector::~SensitiveDetector ( )
override

Definition at line 39 of file SensitiveDetector.cc.

39 {}

Member Function Documentation

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

Definition at line 45 of file SensitiveDetector.cc.

Referenced by SensitiveDetector().

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

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

92  {
93  TrackInformation* info = (TrackInformation*)(aTrack->GetUserInformation());
94  if (!info) {
95  edm::LogWarning("SensitiveDetector") << " no TrackInformation available for trackID= " << aTrack->GetTrackID();
96  throw SimG4Exception("SimG4CoreSensitiveDetector: cannot handle hits for " + GetName());
97  }
98  return info;
99 }
static const TGPicture * info(bool iBackgroundIsBlack)
Log< level::Warning, false > LogWarning
Local3DPoint SensitiveDetector::ConvertToLocal3DPoint ( const G4ThreeVector &  point) const
inlineprotected

Definition at line 46 of file SensitiveDetector.h.

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
void SensitiveDetector::EndOfEvent ( G4HCofThisEvent *  eventHC)
override

Definition at line 43 of file SensitiveDetector.cc.

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

Definition at line 66 of file SensitiveDetector.cc.

References ConvertToLocal3DPoint(), and WorldCoordinates.

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

66  {
67  const G4StepPoint* postStepPoint = step->GetPostStepPoint();
68  const G4ThreeVector& globalCoordinates = postStepPoint->GetPosition();
69  if (cd == WorldCoordinates) {
70  return ConvertToLocal3DPoint(globalCoordinates);
71  }
72  const G4StepPoint* preStepPoint = step->GetPreStepPoint();
73  const G4ThreeVector localCoordinates =
74  preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(globalCoordinates);
75  return ConvertToLocal3DPoint(localCoordinates);
76 }
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
step
Definition: StallMonitor.cc:94
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
void SensitiveDetector::Initialize ( G4HCofThisEvent *  eventHC)
override

Definition at line 41 of file SensitiveDetector.cc.

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

Definition at line 54 of file SensitiveDetector.cc.

References ConvertToLocal3DPoint(), and WorldCoordinates.

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

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

Definition at line 34 of file SensitiveDetector.h.

References m_isCalo.

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

Definition at line 85 of file SensitiveDetector.cc.

References ConvertToLocal3DPoint().

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

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

Definition at line 78 of file SensitiveDetector.cc.

References ConvertToLocal3DPoint().

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

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

Definition at line 106 of file SensitiveDetector.cc.

References edm::isNotFinite().

Referenced by CaloSD::ProcessHits().

106  {
107  const G4Track* CurrentTrk = aStep->GetTrack();
108  const G4ThreeVector& CurrentPos = CurrentTrk->GetPosition();
109  double xyz = CurrentPos.x() + CurrentPos.y() + CurrentPos.z();
110  const G4ThreeVector& CurrentMom = CurrentTrk->GetMomentum();
111  xyz += CurrentMom.x() + CurrentMom.y() + CurrentMom.z();
112 
113  if (edm::isNotFinite(xyz)) {
114  const G4VPhysicalVolume* pCurrentVol = aStep->GetPreStepPoint()->GetPhysicalVolume();
115  const G4String& NameOfVol = pCurrentVol->GetName();
116  throw SimG4Exception("SimG4CoreSensitiveDetector: Corrupted Event - NaN detected in volume " + NameOfVol);
117  }
118 }
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 101 of file SensitiveDetector.cc.

References m_namesOfSD.

Referenced by TkAccumulatingSensitiveDetector::TkAccumulatingSensitiveDetector().

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

Member Data Documentation

bool SensitiveDetector::m_isCalo
private

Definition at line 59 of file SensitiveDetector.h.

Referenced by isCaloSD().

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

Definition at line 58 of file SensitiveDetector.h.

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