CMS 3D CMS Logo

HFWedgeSD.cc
Go to the documentation of this file.
3 
4 #include "G4VPhysicalVolume.hh"
5 #include "G4PVPlacement.hh"
6 #include "G4HCofThisEvent.hh"
7 #include "G4TouchableHistory.hh"
8 #include "G4Track.hh"
9 #include "G4Step.hh"
10 #include "G4VSolid.hh"
11 #include "G4DynamicParticle.hh"
12 #include "G4ParticleDefinition.hh"
13 #include "G4SDManager.hh"
14 #include "G4ios.hh"
15 
16 #include "G4PhysicalConstants.hh"
17 #include "G4SystemOfUnits.hh"
18 
19 HFWedgeSD::HFWedgeSD(const std::string& iname, const DDCompactView & cpv,
20  const SensitiveDetectorCatalog & clg, edm::ParameterSet const & p,
21  const SimTrackManager* manager) :
22  SensitiveCaloDetector(iname, cpv, clg, p),
23  m_trackManager(manager), hcID(-1), theHC(nullptr), currentHit(nullptr) {
24 
25 }
26 
28  delete theHC;
29 }
30 
31 void HFWedgeSD::Initialize(G4HCofThisEvent * HCE) {
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 }
41 
42 G4bool HFWedgeSD::ProcessHits(G4Step * aStep, G4TouchableHistory*) {
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 }
61 
62 void HFWedgeSD::EndOfEvent(G4HCofThisEvent * HCE) {
63 
64  LogDebug("FiberSim") << "HFWedgeSD: Sees" << theHC->entries() << " hits";
65  clear();
66 }
67 
69 
71 
73 
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 }
90 
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 }
112 
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 }
122 
124  hitMap.erase (hitMap.begin(), hitMap.end());
125  previousID = -1;
126 }
127 
128 uint32_t HFWedgeSD::setDetUnitId(const G4Step* aStep) {
129  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
130  return (touch->GetReplicaNumber(0));
131 }
132 
#define LogDebug(id)
HFWedgeSD(const std::string &, const DDCompactView &cpv, const SensitiveDetectorCatalog &clg, edm::ParameterSet const &p, const SimTrackManager *)
Definition: HFWedgeSD.cc:19
std::vector< PCaloHit > PCaloHitContainer
int previousID
Definition: HFWedgeSD.h:52
void updateHit(HFShowerG4Hit *)
Definition: HFWedgeSD.cc:113
void Initialize(G4HCofThisEvent *HCE) override
Definition: HFWedgeSD.cc:31
int hcID
Definition: HFWedgeSD.h:48
G4ThreeVector momDir
Definition: HFWedgeSD.h:54
bool ProcessHits(G4Step *step, G4TouchableHistory *tHistory) override
Definition: HFWedgeSD.cc:42
void setLocalPos(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:46
#define nullptr
void DrawAll() override
Definition: HFWedgeSD.cc:70
Compact representation of the geometrical detector hierarchy.
Definition: DDCompactView.h:80
void fillHits(edm::PCaloHitContainer &, const std::string &) override
Definition: HFWedgeSD.cc:133
const double MeV
std::string const collectionName[nCollections]
Definition: Collections.h:47
double edep
Definition: HFWedgeSD.h:53
void setTrackId(G4int trackId)
Definition: HFShowerG4Hit.h:42
int trackID
Definition: HFWedgeSD.h:52
void clearHits() override
Definition: HFWedgeSD.cc:123
G4bool hitExists()
Definition: HFWedgeSD.cc:74
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 EndOfEvent(G4HCofThisEvent *eventHC) override
Definition: HFWedgeSD.cc:62
void setPrimMomDir(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:48
uint32_t setDetUnitId(const G4Step *) override
Definition: HFWedgeSD.cc:128
HFShowerG4Hit * currentHit
Definition: HFWedgeSD.h:55
void setGlobalPos(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:47
double time
Definition: HFWedgeSD.h:53
void updateEnergy(G4double edep)
Definition: HFShowerG4Hit.h:44
void PrintAll() override
Definition: HFWedgeSD.cc:72
~HFWedgeSD() override
Definition: HFWedgeSD.cc:27
HFShowerG4Hit * createNewHit()
Definition: HFWedgeSD.cc:91
void clear() override
Definition: HFWedgeSD.cc:68
int currentID
Definition: HFWedgeSD.h:52
std::map< int, HFShowerG4Hit * > hitMap
Definition: HFWedgeSD.h:50
G4THitsCollection< HFShowerG4Hit > HFShowerG4HitsCollection
Definition: HFShowerG4Hit.h:59
G4ThreeVector globalPos
Definition: HFWedgeSD.h:54