CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Member Functions | Private Attributes
HFWedgeSD Class Reference

#include <HFWedgeSD.h>

Inheritance diagram for HFWedgeSD:
SensitiveCaloDetector SensitiveDetector

Public Member Functions

void clear () override
 
void clearHits () override
 
void DrawAll () override
 
void EndOfEvent (G4HCofThisEvent *eventHC) override
 
void fillHits (edm::PCaloHitContainer &, const std::string &) override
 
 HFWedgeSD (const std::string &, const SensitiveDetectorCatalog &clg, const SimTrackManager *)
 
void Initialize (G4HCofThisEvent *HCE) override
 
void PrintAll () override
 
bool ProcessHits (G4Step *step, G4TouchableHistory *tHistory) override
 
uint32_t setDetUnitId (const G4Step *) override
 
 ~HFWedgeSD () override
 
- Public Member Functions inherited from SensitiveCaloDetector
virtual void reset ()
 
 SensitiveCaloDetector (const std::string &iname, const SensitiveDetectorCatalog &clg, const std::string &newcollname="")
 
- Public Member Functions inherited from SensitiveDetector
void EndOfEvent (G4HCofThisEvent *eventHC) override
 
const std::vector< std::string > & getNames () const
 
void Initialize (G4HCofThisEvent *eventHC) override
 
bool isCaloSD () const
 
 SensitiveDetector (const std::string &iname, const SensitiveDetectorCatalog &, bool calo, const std::string &newcollname="")
 
 ~SensitiveDetector () override
 

Protected Member Functions

HFShowerG4HitcreateNewHit ()
 
G4bool hitExists ()
 
void updateHit (HFShowerG4Hit *)
 
- Protected Member Functions inherited from SensitiveDetector
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 Attributes

HFShowerG4HitcurrentHit
 
int currentID
 
double edep
 
G4ThreeVector globalPos
 
int hcID
 
std::map< int, HFShowerG4Hit * > hitMap
 
G4ThreeVector localPos
 
G4ThreeVector momDir
 
int previousID
 
HFShowerG4HitsCollectiontheHC
 
double time
 
int trackID
 

Additional Inherited Members

- Protected Types inherited from SensitiveDetector
enum  coordinates { WorldCoordinates, LocalCoordinates }
 

Detailed Description

Definition at line 20 of file HFWedgeSD.h.

Constructor & Destructor Documentation

◆ HFWedgeSD()

HFWedgeSD::HFWedgeSD ( const std::string &  iname,
const SensitiveDetectorCatalog clg,
const SimTrackManager manager 
)
explicit

Definition at line 23 of file HFWedgeSD.cc.

24  : SensitiveCaloDetector(iname, clg), hcID(-1), theHC(nullptr), currentHit(nullptr) {
25  edm::LogVerbatim("FiberSim") << "HFWedgeSD : Instantiated for " << iname;
26 }
Log< level::Info, true > LogVerbatim
int hcID
Definition: HFWedgeSD.h:42
SensitiveCaloDetector(const std::string &iname, const SensitiveDetectorCatalog &clg, const std::string &newcollname="")
HFShowerG4HitsCollection * theHC
Definition: HFWedgeSD.h:43
HFShowerG4Hit * currentHit
Definition: HFWedgeSD.h:49

◆ ~HFWedgeSD()

HFWedgeSD::~HFWedgeSD ( )
override

Definition at line 28 of file HFWedgeSD.cc.

References theHC.

28 { delete theHC; }
HFShowerG4HitsCollection * theHC
Definition: HFWedgeSD.h:43

Member Function Documentation

◆ clear()

void HFWedgeSD::clear ( void  )
override

Definition at line 66 of file HFWedgeSD.cc.

Referenced by EndOfEvent().

66 {}

◆ clearHits()

void HFWedgeSD::clearHits ( )
overridevirtual

Implements SensitiveDetector.

Definition at line 120 of file HFWedgeSD.cc.

References hitMap, and previousID.

Referenced by Initialize().

120  {
121  hitMap.erase(hitMap.begin(), hitMap.end());
122  previousID = -1;
123 }
int previousID
Definition: HFWedgeSD.h:46
std::map< int, HFShowerG4Hit * > hitMap
Definition: HFWedgeSD.h:44

◆ createNewHit()

HFShowerG4Hit * HFWedgeSD::createNewHit ( )
protected

Definition at line 88 of file HFWedgeSD.cc.

References currentID, edep, globalPos, hitMap, localPos, momDir, previousID, HFShowerG4Hit::setGlobalPos(), HFShowerG4Hit::setHitId(), HFShowerG4Hit::setLocalPos(), HFShowerG4Hit::setPrimMomDir(), HFShowerG4Hit::setTime(), HFShowerG4Hit::setTrackId(), theHC, time, trackID, and updateHit().

Referenced by ProcessHits().

88  {
89 #ifdef EDM_ML_DEBUG
90  edm::LogVerbatim("FiberSim") << "HFWedgeSD::CreateNewHit for ID " << currentID << " Track " << trackID
91  << " Edep: " << edep / CLHEP::MeV << " MeV; Time: " << time << " ns; Position (local) "
92  << localPos << " (global ) " << globalPos << " direction " << momDir;
93 #endif
94  HFShowerG4Hit* aHit = new HFShowerG4Hit;
95  aHit->setHitId(currentID);
96  aHit->setTrackId(trackID);
97  aHit->setTime(time);
98  aHit->setLocalPos(localPos);
99  aHit->setGlobalPos(globalPos);
100  aHit->setPrimMomDir(momDir);
101  updateHit(aHit);
102 
103  theHC->insert(aHit);
104  hitMap.insert(std::pair<int, HFShowerG4Hit*>(previousID, aHit));
105 
106  return aHit;
107 }
Log< level::Info, true > LogVerbatim
int previousID
Definition: HFWedgeSD.h:46
void updateHit(HFShowerG4Hit *)
Definition: HFWedgeSD.cc:109
G4ThreeVector momDir
Definition: HFWedgeSD.h:48
void setLocalPos(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:42
std::map< int, HFShowerG4Hit * > hitMap
Definition: HFWedgeSD.h:44
double edep
Definition: HFWedgeSD.h:47
void setTrackId(G4int trackId)
Definition: HFShowerG4Hit.h:38
int trackID
Definition: HFWedgeSD.h:46
G4ThreeVector localPos
Definition: HFWedgeSD.h:48
void setHitId(G4int hitId)
Definition: HFShowerG4Hit.h:37
void setTime(G4double t)
Definition: HFShowerG4Hit.h:41
HFShowerG4HitsCollection * theHC
Definition: HFWedgeSD.h:43
void setPrimMomDir(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:44
void setGlobalPos(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:43
double time
Definition: HFWedgeSD.h:47
int currentID
Definition: HFWedgeSD.h:46
G4ThreeVector globalPos
Definition: HFWedgeSD.h:48

◆ DrawAll()

void HFWedgeSD::DrawAll ( )
override

Definition at line 68 of file HFWedgeSD.cc.

68 {}

◆ EndOfEvent()

void HFWedgeSD::EndOfEvent ( G4HCofThisEvent *  eventHC)
override

Definition at line 61 of file HFWedgeSD.cc.

References clear(), and theHC.

61  {
62  edm::LogVerbatim("FiberSim") << "HFWedgeSD: Sees" << theHC->entries() << " hits";
63  clear();
64 }
Log< level::Info, true > LogVerbatim
HFShowerG4HitsCollection * theHC
Definition: HFWedgeSD.h:43
void clear() override
Definition: HFWedgeSD.cc:66

◆ fillHits()

void HFWedgeSD::fillHits ( edm::PCaloHitContainer ,
const std::string &   
)
overridevirtual

Implements SensitiveCaloDetector.

Definition at line 130 of file HFWedgeSD.cc.

130 {}

◆ hitExists()

G4bool HFWedgeSD::hitExists ( )
protected

Definition at line 72 of file HFWedgeSD.cc.

References currentHit, currentID, hitMap, ALPAKA_ACCELERATOR_NAMESPACE::vertexFinder::it, previousID, and updateHit().

Referenced by ProcessHits().

72  {
73  // Update if in the same detector, time-slice and for same track
74  if (currentID == previousID) {
76  return true;
77  }
78 
79  std::map<int, HFShowerG4Hit*>::const_iterator it = hitMap.find(currentID);
80  if (it != hitMap.end()) {
82  return true;
83  }
84 
85  return false;
86 }
int previousID
Definition: HFWedgeSD.h:46
void updateHit(HFShowerG4Hit *)
Definition: HFWedgeSD.cc:109
std::map< int, HFShowerG4Hit * > hitMap
Definition: HFWedgeSD.h:44
HFShowerG4Hit * currentHit
Definition: HFWedgeSD.h:49
int currentID
Definition: HFWedgeSD.h:46

◆ Initialize()

void HFWedgeSD::Initialize ( G4HCofThisEvent *  HCE)
override

Definition at line 30 of file HFWedgeSD.cc.

References clearHits(), bysipixelclustmulteventfilter_cfi::collectionName, hcID, and theHC.

30  {
31  edm::LogVerbatim("FiberSim") << "HFWedgeSD : Initialize called for " << GetName() << " in collection " << HCE;
32  theHC = new HFShowerG4HitsCollection(GetName(), collectionName[0]);
33  if (hcID < 0)
34  hcID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
35  HCE->AddHitsCollection(hcID, theHC);
36  edm::LogVerbatim("FiberSim") << "HFWedgeSD : Add hit collectrion for " << collectionName[0] << ":" << hcID << ":"
37  << theHC;
38 
39  clearHits();
40 }
Log< level::Info, true > LogVerbatim
int hcID
Definition: HFWedgeSD.h:42
void clearHits() override
Definition: HFWedgeSD.cc:120
HFShowerG4HitsCollection * theHC
Definition: HFWedgeSD.h:43
G4THitsCollection< HFShowerG4Hit > HFShowerG4HitsCollection
Definition: HFShowerG4Hit.h:55

◆ PrintAll()

void HFWedgeSD::PrintAll ( )
override

Definition at line 70 of file HFWedgeSD.cc.

70 {}

◆ ProcessHits()

G4bool HFWedgeSD::ProcessHits ( G4Step *  step,
G4TouchableHistory *  tHistory 
)
overridevirtual

Implements SensitiveDetector.

Definition at line 42 of file HFWedgeSD.cc.

References createNewHit(), currentHit, currentID, edep, globalPos, hitExists(), localPos, momDir, setDetUnitId(), time, and trackID.

42  {
43  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
44  const G4VTouchable* touch = preStepPoint->GetTouchable();
45  currentID = setDetUnitId(aStep);
46  trackID = aStep->GetTrack()->GetTrackID();
47  edep = aStep->GetTotalEnergyDeposit();
48  time = (preStepPoint->GetGlobalTime()) / ns;
49 
50  globalPos = preStepPoint->GetPosition();
51  localPos = touch->GetHistory()->GetTopTransform().TransformPoint(globalPos);
52  const G4DynamicParticle* particle = aStep->GetTrack()->GetDynamicParticle();
53  momDir = particle->GetMomentumDirection();
54 
55  if (hitExists() == false && edep > 0.)
57 
58  return true;
59 }
G4ThreeVector momDir
Definition: HFWedgeSD.h:48
double edep
Definition: HFWedgeSD.h:47
int trackID
Definition: HFWedgeSD.h:46
G4bool hitExists()
Definition: HFWedgeSD.cc:72
G4ThreeVector localPos
Definition: HFWedgeSD.h:48
uint32_t setDetUnitId(const G4Step *) override
Definition: HFWedgeSD.cc:125
HFShowerG4Hit * currentHit
Definition: HFWedgeSD.h:49
double time
Definition: HFWedgeSD.h:47
HFShowerG4Hit * createNewHit()
Definition: HFWedgeSD.cc:88
int currentID
Definition: HFWedgeSD.h:46
G4ThreeVector globalPos
Definition: HFWedgeSD.h:48

◆ setDetUnitId()

uint32_t HFWedgeSD::setDetUnitId ( const G4Step *  aStep)
overridevirtual

Implements SensitiveDetector.

Definition at line 125 of file HFWedgeSD.cc.

Referenced by ProcessHits().

125  {
126  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
127  return (touch->GetReplicaNumber(0));
128 }

◆ updateHit()

void HFWedgeSD::updateHit ( HFShowerG4Hit aHit)
protected

Definition at line 109 of file HFWedgeSD.cc.

References currentID, edep, previousID, and HFShowerG4Hit::updateEnergy().

Referenced by createNewHit(), and hitExists().

109  {
110  if (edep != 0) {
111  aHit->updateEnergy(edep);
112 #ifdef EDM_ML_DEBUG
113  edm::LogVerbatim("FiberSim") << "HFWedgeSD: Add energy deposit in " << currentID << " edep " << edep / CLHEP::MeV
114  << " MeV";
115 #endif
116  }
118 }
Log< level::Info, true > LogVerbatim
int previousID
Definition: HFWedgeSD.h:46
double edep
Definition: HFWedgeSD.h:47
void updateEnergy(G4double edep)
Definition: HFShowerG4Hit.h:40
int currentID
Definition: HFWedgeSD.h:46

Member Data Documentation

◆ currentHit

HFShowerG4Hit* HFWedgeSD::currentHit
private

Definition at line 49 of file HFWedgeSD.h.

Referenced by hitExists(), and ProcessHits().

◆ currentID

int HFWedgeSD::currentID
private

Definition at line 46 of file HFWedgeSD.h.

Referenced by createNewHit(), hitExists(), ProcessHits(), and updateHit().

◆ edep

double HFWedgeSD::edep
private

Definition at line 47 of file HFWedgeSD.h.

Referenced by createNewHit(), ProcessHits(), and updateHit().

◆ globalPos

G4ThreeVector HFWedgeSD::globalPos
private

Definition at line 48 of file HFWedgeSD.h.

Referenced by createNewHit(), and ProcessHits().

◆ hcID

int HFWedgeSD::hcID
private

Definition at line 42 of file HFWedgeSD.h.

Referenced by Initialize().

◆ hitMap

std::map<int, HFShowerG4Hit*> HFWedgeSD::hitMap
private

Definition at line 44 of file HFWedgeSD.h.

Referenced by clearHits(), createNewHit(), and hitExists().

◆ localPos

G4ThreeVector HFWedgeSD::localPos
private

Definition at line 48 of file HFWedgeSD.h.

Referenced by createNewHit(), and ProcessHits().

◆ momDir

G4ThreeVector HFWedgeSD::momDir
private

Definition at line 48 of file HFWedgeSD.h.

Referenced by createNewHit(), and ProcessHits().

◆ previousID

int HFWedgeSD::previousID
private

Definition at line 46 of file HFWedgeSD.h.

Referenced by clearHits(), createNewHit(), hitExists(), and updateHit().

◆ theHC

HFShowerG4HitsCollection* HFWedgeSD::theHC
private

Definition at line 43 of file HFWedgeSD.h.

Referenced by createNewHit(), EndOfEvent(), Initialize(), and ~HFWedgeSD().

◆ time

double HFWedgeSD::time
private

Definition at line 47 of file HFWedgeSD.h.

Referenced by createNewHit(), and ProcessHits().

◆ trackID

int HFWedgeSD::trackID
private

Definition at line 46 of file HFWedgeSD.h.

Referenced by createNewHit(), and ProcessHits().