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 DDCompactView &cpv, const SensitiveDetectorCatalog &clg, edm::ParameterSet const &p, 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 DDCompactView &cpv, const SensitiveDetectorCatalog &clg, edm::ParameterSet const &p)
 
- 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 DDCompactView &cpv, const SensitiveDetectorCatalog &, edm::ParameterSet const &p, bool calo)
 
 ~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
 
const SimTrackManagerm_trackManager
 
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 18 of file HFWedgeSD.h.

Constructor & Destructor Documentation

HFWedgeSD::HFWedgeSD ( const std::string &  iname,
const DDCompactView cpv,
const SensitiveDetectorCatalog clg,
edm::ParameterSet const &  p,
const SimTrackManager manager 
)
explicit

Definition at line 19 of file HFWedgeSD.cc.

21  :
22  SensitiveCaloDetector(iname, cpv, clg, p),
23  m_trackManager(manager), hcID(-1), theHC(nullptr), currentHit(nullptr) {
24 
25 }
SensitiveCaloDetector(const std::string &iname, const DDCompactView &cpv, const SensitiveDetectorCatalog &clg, edm::ParameterSet const &p)
int hcID
Definition: HFWedgeSD.h:48
const SimTrackManager * m_trackManager
Definition: HFWedgeSD.h:46
HFShowerG4HitsCollection * theHC
Definition: HFWedgeSD.h:49
HFShowerG4Hit * currentHit
Definition: HFWedgeSD.h:55
HFWedgeSD::~HFWedgeSD ( )
override

Definition at line 27 of file HFWedgeSD.cc.

References theHC.

27  {
28  delete theHC;
29 }
HFShowerG4HitsCollection * theHC
Definition: HFWedgeSD.h:49

Member Function Documentation

void HFWedgeSD::clear ( void  )
override

Definition at line 68 of file HFWedgeSD.cc.

Referenced by EndOfEvent().

68 {}
void HFWedgeSD::clearHits ( )
overridevirtual

Implements SensitiveDetector.

Definition at line 123 of file HFWedgeSD.cc.

References hitMap, and previousID.

Referenced by Initialize().

123  {
124  hitMap.erase (hitMap.begin(), hitMap.end());
125  previousID = -1;
126 }
int previousID
Definition: HFWedgeSD.h:52
std::map< int, HFShowerG4Hit * > hitMap
Definition: HFWedgeSD.h:50
HFShowerG4Hit * HFWedgeSD::createNewHit ( )
protected

Definition at line 91 of file HFWedgeSD.cc.

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

Referenced by ProcessHits().

91  {
92 
93  LogDebug("FiberSim") << "HFWedgeSD::CreateNewHit for ID " << currentID
94  << " Track " << trackID << " Edep: " << edep/MeV
95  << " MeV; Time: " << time << " ns; Position (local) "
96  << localPos << " (global ) " << globalPos
97  << " direction " << momDir;
98  HFShowerG4Hit* aHit = new HFShowerG4Hit;
99  aHit->setHitId(currentID);
100  aHit->setTrackId(trackID);
101  aHit->setTime(time);
102  aHit->setLocalPos(localPos);
103  aHit->setGlobalPos(globalPos);
104  aHit->setPrimMomDir(momDir);
105  updateHit(aHit);
106 
107  theHC->insert(aHit);
108  hitMap.insert(std::pair<int,HFShowerG4Hit*>(previousID,aHit));
109 
110  return aHit;
111 }
#define LogDebug(id)
int previousID
Definition: HFWedgeSD.h:52
void updateHit(HFShowerG4Hit *)
Definition: HFWedgeSD.cc:113
G4ThreeVector momDir
Definition: HFWedgeSD.h:54
void setLocalPos(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:46
const double MeV
double edep
Definition: HFWedgeSD.h:53
void setTrackId(G4int trackId)
Definition: HFShowerG4Hit.h:42
int trackID
Definition: HFWedgeSD.h:52
G4ThreeVector localPos
Definition: HFWedgeSD.h:54
void setHitId(G4int hitId)
Definition: HFShowerG4Hit.h:41
void setTime(G4double t)
Definition: HFShowerG4Hit.h:45
HFShowerG4HitsCollection * theHC
Definition: HFWedgeSD.h:49
void setPrimMomDir(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:48
void setGlobalPos(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:47
double time
Definition: HFWedgeSD.h:53
int currentID
Definition: HFWedgeSD.h:52
std::map< int, HFShowerG4Hit * > hitMap
Definition: HFWedgeSD.h:50
G4ThreeVector globalPos
Definition: HFWedgeSD.h:54
void HFWedgeSD::DrawAll ( )
override

Definition at line 70 of file HFWedgeSD.cc.

70 {}
void HFWedgeSD::EndOfEvent ( G4HCofThisEvent *  eventHC)
override

Definition at line 62 of file HFWedgeSD.cc.

References clear(), LogDebug, and theHC.

62  {
63 
64  LogDebug("FiberSim") << "HFWedgeSD: Sees" << theHC->entries() << " hits";
65  clear();
66 }
#define LogDebug(id)
HFShowerG4HitsCollection * theHC
Definition: HFWedgeSD.h:49
void clear() override
Definition: HFWedgeSD.cc:68
void HFWedgeSD::fillHits ( edm::PCaloHitContainer ,
const std::string &   
)
overridevirtual

Implements SensitiveCaloDetector.

Definition at line 133 of file HFWedgeSD.cc.

133 {}
G4bool HFWedgeSD::hitExists ( )
protected

Definition at line 74 of file HFWedgeSD.cc.

References currentHit, currentID, hitMap, previousID, and updateHit().

Referenced by ProcessHits().

74  {
75 
76  // Update if in the same detector, time-slice and for same track
77  if (currentID == previousID) {
79  return true;
80  }
81 
82  std::map<int,HFShowerG4Hit*>::const_iterator it = hitMap.find(currentID);
83  if (it != hitMap.end()) {
85  return true;
86  }
87 
88  return false;
89 }
int previousID
Definition: HFWedgeSD.h:52
void updateHit(HFShowerG4Hit *)
Definition: HFWedgeSD.cc:113
HFShowerG4Hit * currentHit
Definition: HFWedgeSD.h:55
int currentID
Definition: HFWedgeSD.h:52
std::map< int, HFShowerG4Hit * > hitMap
Definition: HFWedgeSD.h:50
void HFWedgeSD::Initialize ( G4HCofThisEvent *  HCE)
override

Definition at line 31 of file HFWedgeSD.cc.

References clearHits(), ecaldqm::collectionName, hcID, LogDebug, and theHC.

31  {
32 
33  LogDebug("FiberSim") << "HFWedgeSD : Initialize called for " << GetName();
34  theHC = new HFShowerG4HitsCollection(GetName(), collectionName[0]);
35  if (hcID<0)
36  hcID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
37  HCE->AddHitsCollection(hcID, theHC);
38 
39  clearHits();
40 }
#define LogDebug(id)
int hcID
Definition: HFWedgeSD.h:48
std::string const collectionName[nCollections]
Definition: Collections.h:47
void clearHits() override
Definition: HFWedgeSD.cc:123
HFShowerG4HitsCollection * theHC
Definition: HFWedgeSD.h:49
G4THitsCollection< HFShowerG4Hit > HFShowerG4HitsCollection
Definition: HFShowerG4Hit.h:59
void HFWedgeSD::PrintAll ( )
override

Definition at line 72 of file HFWedgeSD.cc.

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

Implements SensitiveDetector.

Definition at line 128 of file HFWedgeSD.cc.

Referenced by ProcessHits().

128  {
129  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
130  return (touch->GetReplicaNumber(0));
131 }
void HFWedgeSD::updateHit ( HFShowerG4Hit aHit)
protected

Definition at line 113 of file HFWedgeSD.cc.

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

Referenced by createNewHit(), and hitExists().

113  {
114 
115  if (edep != 0) {
116  aHit->updateEnergy(edep);
117  LogDebug("FiberSim") << "HFWedgeSD: Add energy deposit in " << currentID
118  << " edep " << edep/MeV << " MeV";
119  }
121 }
#define LogDebug(id)
int previousID
Definition: HFWedgeSD.h:52
const double MeV
double edep
Definition: HFWedgeSD.h:53
void updateEnergy(G4double edep)
Definition: HFShowerG4Hit.h:44
int currentID
Definition: HFWedgeSD.h:52

Member Data Documentation

HFShowerG4Hit* HFWedgeSD::currentHit
private

Definition at line 55 of file HFWedgeSD.h.

Referenced by hitExists(), and ProcessHits().

int HFWedgeSD::currentID
private

Definition at line 52 of file HFWedgeSD.h.

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

double HFWedgeSD::edep
private

Definition at line 53 of file HFWedgeSD.h.

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

G4ThreeVector HFWedgeSD::globalPos
private

Definition at line 54 of file HFWedgeSD.h.

Referenced by createNewHit(), and ProcessHits().

int HFWedgeSD::hcID
private

Definition at line 48 of file HFWedgeSD.h.

Referenced by Initialize().

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

Definition at line 50 of file HFWedgeSD.h.

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

G4ThreeVector HFWedgeSD::localPos
private

Definition at line 54 of file HFWedgeSD.h.

Referenced by createNewHit(), and ProcessHits().

const SimTrackManager* HFWedgeSD::m_trackManager
private

Definition at line 46 of file HFWedgeSD.h.

G4ThreeVector HFWedgeSD::momDir
private

Definition at line 54 of file HFWedgeSD.h.

Referenced by createNewHit(), and ProcessHits().

int HFWedgeSD::previousID
private

Definition at line 52 of file HFWedgeSD.h.

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

HFShowerG4HitsCollection* HFWedgeSD::theHC
private

Definition at line 49 of file HFWedgeSD.h.

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

double HFWedgeSD::time
private

Definition at line 53 of file HFWedgeSD.h.

Referenced by createNewHit(), and ProcessHits().

int HFWedgeSD::trackID
private

Definition at line 52 of file HFWedgeSD.h.

Referenced by createNewHit(), and ProcessHits().