CMS 3D CMS Logo

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

virtual void clear ()
 
virtual void DrawAll ()
 
virtual void EndOfEvent (G4HCofThisEvent *eventHC)
 
 HFWedgeSD (std::string name, const DDCompactView &cpv, const SensitiveDetectorCatalog &clg, edm::ParameterSet const &p, const SimTrackManager *)
 
virtual void Initialize (G4HCofThisEvent *HCE)
 
virtual void PrintAll ()
 
virtual bool ProcessHits (G4Step *step, G4TouchableHistory *tHistory)
 
virtual ~HFWedgeSD ()
 
- Public Member Functions inherited from SensitiveCaloDetector
 SensitiveCaloDetector (std::string &iname, const DDCompactView &cpv, const SensitiveDetectorCatalog &clg, edm::ParameterSet const &p)
 
- Public Member Functions inherited from SensitiveDetector
virtual void AssignSD (const std::string &vname)
 
Local3DPoint ConvertToLocal3DPoint (const G4ThreeVector &point)
 
Local3DPoint FinalStepPosition (G4Step *s, coordinates)
 
virtual std::vector< std::string > getNames ()
 
Local3DPoint InitialStepPosition (G4Step *s, coordinates)
 
std::string nameOfSD ()
 
void NaNTrap (G4Step *step)
 
void Register ()
 
 SensitiveDetector (std::string &iname, const DDCompactView &cpv, const SensitiveDetectorCatalog &, edm::ParameterSet const &p)
 
virtual ~SensitiveDetector ()
 

Protected Member Functions

virtual void clearHits ()
 
HFShowerG4HitcreateNewHit ()
 
virtual void fillHits (edm::PCaloHitContainer &, std::string)
 
G4bool hitExists ()
 
virtual uint32_t setDetUnitId (G4Step *)
 
void updateHit (HFShowerG4Hit *)
 

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
 
std::string theName
 
double time
 
int trackID
 

Additional Inherited Members

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

Detailed Description

Definition at line 22 of file HFWedgeSD.h.

Constructor & Destructor Documentation

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

Definition at line 19 of file HFWedgeSD.cc.

References SensitiveDetector::AssignSD(), ecaldqm::collectionName, LogDebug, SensitiveDetectorCatalog::logicalNames(), and SensitiveDetector::Register().

21  :
22  SensitiveCaloDetector(name, cpv, clg, p), theName(name),
23  m_trackManager(manager), hcID(-1), theHC(0), currentHit(0) {
24 
25  collectionName.insert(name);
26  LogDebug("FiberSim") << "***************************************************"
27  << "\n"
28  << "* *"
29  << "\n"
30  << "* Constructing a HFWedgeSD with name " << GetName()
31  << "\n"
32  << "* *"
33  << "\n"
34  << "***************************************************";
35  //
36  // Now attach the right detectors (LogicalVolumes) to me
37  //
38  const std::vector<std::string>& lvNames = clg.logicalNames(name);
39  this->Register();
40  for (std::vector<std::string>::const_iterator it=lvNames.begin();
41  it !=lvNames.end(); it++){
42  this->AssignSD(*it);
43  LogDebug("FiberSim") << "HFWedgeSD : Assigns SD to LV " << (*it);
44  }
45 }
#define LogDebug(id)
int hcID
Definition: HFWedgeSD.h:54
SensitiveCaloDetector(std::string &iname, const DDCompactView &cpv, const SensitiveDetectorCatalog &clg, edm::ParameterSet const &p)
const std::vector< std::string > & logicalNames(const std::string &readoutName) const
virtual void AssignSD(const std::string &vname)
std::string const collectionName[nCollections]
Definition: Collections.h:45
const SimTrackManager * m_trackManager
Definition: HFWedgeSD.h:52
HFShowerG4HitsCollection * theHC
Definition: HFWedgeSD.h:55
HFShowerG4Hit * currentHit
Definition: HFWedgeSD.h:61
std::string theName
Definition: HFWedgeSD.h:51
HFWedgeSD::~HFWedgeSD ( )
virtual

Definition at line 47 of file HFWedgeSD.cc.

References theHC.

47  {
48  if (theHC) delete theHC;
49 }
HFShowerG4HitsCollection * theHC
Definition: HFWedgeSD.h:55

Member Function Documentation

void HFWedgeSD::clear ( void  )
virtual

Definition at line 88 of file HFWedgeSD.cc.

Referenced by EndOfEvent().

88 {}
void HFWedgeSD::clearHits ( )
protectedvirtual

Implements SensitiveDetector.

Definition at line 143 of file HFWedgeSD.cc.

References hitMap, and previousID.

Referenced by Initialize().

143  {
144  hitMap.erase (hitMap.begin(), hitMap.end());
145  previousID = -1;
146 }
int previousID
Definition: HFWedgeSD.h:58
std::map< int, HFShowerG4Hit * > hitMap
Definition: HFWedgeSD.h:56
HFShowerG4Hit * HFWedgeSD::createNewHit ( )
protected

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

111  {
112 
113  LogDebug("FiberSim") << "HFWedgeSD::CreateNewHit for ID " << currentID
114  << " Track " << trackID << " Edep: " << edep/MeV
115  << " MeV; Time: " << time << " ns; Position (local) "
116  << localPos << " (global ) " << globalPos
117  << " direction " << momDir;
118  HFShowerG4Hit* aHit = new HFShowerG4Hit;
119  aHit->setHitId(currentID);
120  aHit->setTrackId(trackID);
121  aHit->setTime(time);
122  aHit->setLocalPos(localPos);
123  aHit->setGlobalPos(globalPos);
124  aHit->setPrimMomDir(momDir);
125  updateHit(aHit);
126 
127  theHC->insert(aHit);
128  hitMap.insert(std::pair<int,HFShowerG4Hit*>(previousID,aHit));
129 
130  return aHit;
131 }
#define LogDebug(id)
int previousID
Definition: HFWedgeSD.h:58
void updateHit(HFShowerG4Hit *)
Definition: HFWedgeSD.cc:133
G4ThreeVector momDir
Definition: HFWedgeSD.h:60
void setLocalPos(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:47
const double MeV
double edep
Definition: HFWedgeSD.h:59
void setTrackId(G4int trackId)
Definition: HFShowerG4Hit.h:43
int trackID
Definition: HFWedgeSD.h:58
G4ThreeVector localPos
Definition: HFWedgeSD.h:60
void setHitId(G4int hitId)
Definition: HFShowerG4Hit.h:42
void setTime(G4double t)
Definition: HFShowerG4Hit.h:46
HFShowerG4HitsCollection * theHC
Definition: HFWedgeSD.h:55
void setPrimMomDir(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:49
void setGlobalPos(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:48
double time
Definition: HFWedgeSD.h:59
int currentID
Definition: HFWedgeSD.h:58
std::map< int, HFShowerG4Hit * > hitMap
Definition: HFWedgeSD.h:56
G4ThreeVector globalPos
Definition: HFWedgeSD.h:60
void HFWedgeSD::DrawAll ( )
virtual

Definition at line 90 of file HFWedgeSD.cc.

90 {}
void HFWedgeSD::EndOfEvent ( G4HCofThisEvent *  eventHC)
virtual

Reimplemented from SensitiveDetector.

Definition at line 82 of file HFWedgeSD.cc.

References clear(), LogDebug, and theHC.

82  {
83 
84  LogDebug("FiberSim") << "HFWedgeSD: Sees" << theHC->entries() << " hits";
85  clear();
86 }
#define LogDebug(id)
virtual void clear()
Definition: HFWedgeSD.cc:88
HFShowerG4HitsCollection * theHC
Definition: HFWedgeSD.h:55
void HFWedgeSD::fillHits ( edm::PCaloHitContainer ,
std::string   
)
protectedvirtual

Implements SensitiveCaloDetector.

Definition at line 153 of file HFWedgeSD.cc.

153 {}
G4bool HFWedgeSD::hitExists ( )
protected

Definition at line 94 of file HFWedgeSD.cc.

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

Referenced by ProcessHits().

94  {
95 
96  // Update if in the same detector, time-slice and for same track
97  if (currentID == previousID) {
99  return true;
100  }
101 
102  std::map<int,HFShowerG4Hit*>::const_iterator it = hitMap.find(currentID);
103  if (it != hitMap.end()) {
105  return true;
106  }
107 
108  return false;
109 }
int previousID
Definition: HFWedgeSD.h:58
void updateHit(HFShowerG4Hit *)
Definition: HFWedgeSD.cc:133
HFShowerG4Hit * currentHit
Definition: HFWedgeSD.h:61
int currentID
Definition: HFWedgeSD.h:58
std::map< int, HFShowerG4Hit * > hitMap
Definition: HFWedgeSD.h:56
void HFWedgeSD::Initialize ( G4HCofThisEvent *  HCE)
virtual

Reimplemented from SensitiveDetector.

Definition at line 51 of file HFWedgeSD.cc.

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

51  {
52 
53  LogDebug("FiberSim") << "HFWedgeSD : Initialize called for " << GetName();
54  theHC = new HFShowerG4HitsCollection(GetName(), collectionName[0]);
55  if (hcID<0)
56  hcID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
57  HCE->AddHitsCollection(hcID, theHC);
58 
59  clearHits();
60 }
#define LogDebug(id)
virtual void clearHits()
Definition: HFWedgeSD.cc:143
int hcID
Definition: HFWedgeSD.h:54
std::string const collectionName[nCollections]
Definition: Collections.h:45
HFShowerG4HitsCollection * theHC
Definition: HFWedgeSD.h:55
G4THitsCollection< HFShowerG4Hit > HFShowerG4HitsCollection
Definition: HFShowerG4Hit.h:60
void HFWedgeSD::PrintAll ( )
virtual

Definition at line 92 of file HFWedgeSD.cc.

92 {}
G4bool HFWedgeSD::ProcessHits ( G4Step *  step,
G4TouchableHistory *  tHistory 
)
virtual

Implements SensitiveDetector.

Definition at line 62 of file HFWedgeSD.cc.

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

62  {
63 
64  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
65  const G4VTouchable* touch = preStepPoint->GetTouchable();
66  currentID = setDetUnitId(aStep);
67  trackID = aStep->GetTrack()->GetTrackID();
68  edep = aStep->GetTotalEnergyDeposit();
69  time = (preStepPoint->GetGlobalTime())/ns;
70 
71  globalPos = preStepPoint->GetPosition();
72  localPos = touch->GetHistory()->GetTopTransform().TransformPoint(globalPos);
73  const G4DynamicParticle* particle = aStep->GetTrack()->GetDynamicParticle();
74  momDir = particle->GetMomentumDirection();
75 
76  if (hitExists() == false && edep>0.)
78 
79  return true;
80 }
G4ThreeVector momDir
Definition: HFWedgeSD.h:60
virtual uint32_t setDetUnitId(G4Step *)
Definition: HFWedgeSD.cc:148
double edep
Definition: HFWedgeSD.h:59
int trackID
Definition: HFWedgeSD.h:58
G4bool hitExists()
Definition: HFWedgeSD.cc:94
G4ThreeVector localPos
Definition: HFWedgeSD.h:60
HFShowerG4Hit * currentHit
Definition: HFWedgeSD.h:61
double time
Definition: HFWedgeSD.h:59
HFShowerG4Hit * createNewHit()
Definition: HFWedgeSD.cc:111
int currentID
Definition: HFWedgeSD.h:58
G4ThreeVector globalPos
Definition: HFWedgeSD.h:60
uint32_t HFWedgeSD::setDetUnitId ( G4Step *  aStep)
protectedvirtual

Implements SensitiveDetector.

Definition at line 148 of file HFWedgeSD.cc.

Referenced by ProcessHits().

148  {
149  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
150  return (touch->GetReplicaNumber(0));
151 }
void HFWedgeSD::updateHit ( HFShowerG4Hit aHit)
protected

Definition at line 133 of file HFWedgeSD.cc.

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

Referenced by createNewHit(), and hitExists().

133  {
134 
135  if (edep != 0) {
136  aHit->updateEnergy(edep);
137  LogDebug("FiberSim") << "HFWedgeSD: Add energy deposit in " << currentID
138  << " edep " << edep/MeV << " MeV";
139  }
141 }
#define LogDebug(id)
int previousID
Definition: HFWedgeSD.h:58
const double MeV
double edep
Definition: HFWedgeSD.h:59
void updateEnergy(G4double edep)
Definition: HFShowerG4Hit.h:45
int currentID
Definition: HFWedgeSD.h:58

Member Data Documentation

HFShowerG4Hit* HFWedgeSD::currentHit
private

Definition at line 61 of file HFWedgeSD.h.

Referenced by hitExists(), and ProcessHits().

int HFWedgeSD::currentID
private

Definition at line 58 of file HFWedgeSD.h.

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

double HFWedgeSD::edep
private

Definition at line 59 of file HFWedgeSD.h.

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

G4ThreeVector HFWedgeSD::globalPos
private

Definition at line 60 of file HFWedgeSD.h.

Referenced by createNewHit(), and ProcessHits().

int HFWedgeSD::hcID
private

Definition at line 54 of file HFWedgeSD.h.

Referenced by Initialize().

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

Definition at line 56 of file HFWedgeSD.h.

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

G4ThreeVector HFWedgeSD::localPos
private

Definition at line 60 of file HFWedgeSD.h.

Referenced by createNewHit(), and ProcessHits().

const SimTrackManager* HFWedgeSD::m_trackManager
private

Definition at line 52 of file HFWedgeSD.h.

G4ThreeVector HFWedgeSD::momDir
private

Definition at line 60 of file HFWedgeSD.h.

Referenced by createNewHit(), and ProcessHits().

int HFWedgeSD::previousID
private

Definition at line 58 of file HFWedgeSD.h.

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

HFShowerG4HitsCollection* HFWedgeSD::theHC
private

Definition at line 55 of file HFWedgeSD.h.

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

std::string HFWedgeSD::theName
private

Definition at line 51 of file HFWedgeSD.h.

double HFWedgeSD::time
private

Definition at line 59 of file HFWedgeSD.h.

Referenced by createNewHit(), and ProcessHits().

int HFWedgeSD::trackID
private

Definition at line 58 of file HFWedgeSD.h.

Referenced by createNewHit(), and ProcessHits().