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 //#define EDM_ML_DEBUG
20 
22  : SensitiveCaloDetector(iname, clg), hcID(-1), theHC(nullptr), currentHit(nullptr) {
23  edm::LogVerbatim("FiberSim") << "HFWedgeSD : Instantiated for " << iname;
24 }
25 
27 
28 void HFWedgeSD::Initialize(G4HCofThisEvent* HCE) {
29  edm::LogVerbatim("FiberSim") << "HFWedgeSD : Initialize called for " << GetName() << " in collection " << HCE;
30  theHC = new HFShowerG4HitsCollection(GetName(), collectionName[0]);
31  if (hcID < 0)
32  hcID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
33  HCE->AddHitsCollection(hcID, theHC);
34  edm::LogVerbatim("FiberSim") << "HFWedgeSD : Add hit collectrion for " << collectionName[0] << ":" << hcID << ":"
35  << theHC;
36 
37  clearHits();
38 }
39 
40 G4bool HFWedgeSD::ProcessHits(G4Step* aStep, G4TouchableHistory*) {
41  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
42  const G4VTouchable* touch = preStepPoint->GetTouchable();
43  currentID = setDetUnitId(aStep);
44  trackID = aStep->GetTrack()->GetTrackID();
45  edep = aStep->GetTotalEnergyDeposit();
46  time = (preStepPoint->GetGlobalTime()) / ns;
47 
48  globalPos = preStepPoint->GetPosition();
49  localPos = touch->GetHistory()->GetTopTransform().TransformPoint(globalPos);
50  const G4DynamicParticle* particle = aStep->GetTrack()->GetDynamicParticle();
51  momDir = particle->GetMomentumDirection();
52 
53  if (hitExists() == false && edep > 0.)
55 
56  return true;
57 }
58 
59 void HFWedgeSD::EndOfEvent(G4HCofThisEvent* HCE) {
60  edm::LogVerbatim("FiberSim") << "HFWedgeSD: Sees" << theHC->entries() << " hits";
61  clear();
62 }
63 
65 
67 
69 
71  // Update if in the same detector, time-slice and for same track
72  if (currentID == previousID) {
74  return true;
75  }
76 
77  std::map<int, HFShowerG4Hit*>::const_iterator it = hitMap.find(currentID);
78  if (it != hitMap.end()) {
80  return true;
81  }
82 
83  return false;
84 }
85 
87 #ifdef EDM_ML_DEBUG
88  edm::LogVerbatim("FiberSim") << "HFWedgeSD::CreateNewHit for ID " << currentID << " Track " << trackID
89  << " Edep: " << edep / CLHEP::MeV << " MeV; Time: " << time << " ns; Position (local) "
90  << localPos << " (global ) " << globalPos << " direction " << momDir;
91 #endif
92  HFShowerG4Hit* aHit = new HFShowerG4Hit;
93  aHit->setHitId(currentID);
94  aHit->setTrackId(trackID);
95  aHit->setTime(time);
96  aHit->setLocalPos(localPos);
97  aHit->setGlobalPos(globalPos);
98  aHit->setPrimMomDir(momDir);
99  updateHit(aHit);
100 
101  theHC->insert(aHit);
102  hitMap.insert(std::pair<int, HFShowerG4Hit*>(previousID, aHit));
103 
104  return aHit;
105 }
106 
108  if (edep != 0) {
109  aHit->updateEnergy(edep);
110 #ifdef EDM_ML_DEBUG
111  edm::LogVerbatim("FiberSim") << "HFWedgeSD: Add energy deposit in " << currentID << " edep " << edep / CLHEP::MeV
112  << " MeV";
113 #endif
114  }
116 }
117 
119  hitMap.erase(hitMap.begin(), hitMap.end());
120  previousID = -1;
121 }
122 
123 uint32_t HFWedgeSD::setDetUnitId(const G4Step* aStep) {
124  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
125  return (touch->GetReplicaNumber(0));
126 }
127 
Log< level::Info, true > LogVerbatim
HFWedgeSD(const std::string &, const SensitiveDetectorCatalog &clg, const SimTrackManager *)
Definition: HFWedgeSD.cc:21
std::vector< PCaloHit > PCaloHitContainer
int previousID
Definition: HFWedgeSD.h:46
void updateHit(HFShowerG4Hit *)
Definition: HFWedgeSD.cc:107
void Initialize(G4HCofThisEvent *HCE) override
Definition: HFWedgeSD.cc:28
int hcID
Definition: HFWedgeSD.h:42
G4ThreeVector momDir
Definition: HFWedgeSD.h:48
bool ProcessHits(G4Step *step, G4TouchableHistory *tHistory) override
Definition: HFWedgeSD.cc:40
void setLocalPos(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:42
void DrawAll() override
Definition: HFWedgeSD.cc:66
void fillHits(edm::PCaloHitContainer &, const std::string &) override
Definition: HFWedgeSD.cc:128
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
void clearHits() override
Definition: HFWedgeSD.cc:118
G4bool hitExists()
Definition: HFWedgeSD.cc:70
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 EndOfEvent(G4HCofThisEvent *eventHC) override
Definition: HFWedgeSD.cc:59
void setPrimMomDir(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:44
uint32_t setDetUnitId(const G4Step *) override
Definition: HFWedgeSD.cc:123
HFShowerG4Hit * currentHit
Definition: HFWedgeSD.h:49
void setGlobalPos(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:43
double time
Definition: HFWedgeSD.h:47
void updateEnergy(G4double edep)
Definition: HFShowerG4Hit.h:40
void PrintAll() override
Definition: HFWedgeSD.cc:68
~HFWedgeSD() override
Definition: HFWedgeSD.cc:26
HFShowerG4Hit * createNewHit()
Definition: HFWedgeSD.cc:86
void clear() override
Definition: HFWedgeSD.cc:64
int currentID
Definition: HFWedgeSD.h:46
G4THitsCollection< HFShowerG4Hit > HFShowerG4HitsCollection
Definition: HFShowerG4Hit.h:55
G4ThreeVector globalPos
Definition: HFWedgeSD.h:48