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 
18 #include <CLHEP/Units/SystemOfUnits.h>
19 using CLHEP::ns;
20 
21 //#define EDM_ML_DEBUG
22 
24  : SensitiveCaloDetector(iname, clg), hcID(-1), theHC(nullptr), currentHit(nullptr) {
25  edm::LogVerbatim("FiberSim") << "HFWedgeSD : Instantiated for " << iname;
26 }
27 
29 
30 void HFWedgeSD::Initialize(G4HCofThisEvent* HCE) {
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 }
41 
42 G4bool HFWedgeSD::ProcessHits(G4Step* aStep, G4TouchableHistory*) {
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 }
60 
61 void HFWedgeSD::EndOfEvent(G4HCofThisEvent* HCE) {
62  edm::LogVerbatim("FiberSim") << "HFWedgeSD: Sees" << theHC->entries() << " hits";
63  clear();
64 }
65 
67 
69 
71 
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 }
87 
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 }
108 
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 }
119 
121  hitMap.erase(hitMap.begin(), hitMap.end());
122  previousID = -1;
123 }
124 
125 uint32_t HFWedgeSD::setDetUnitId(const G4Step* aStep) {
126  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
127  return (touch->GetReplicaNumber(0));
128 }
129 
Log< level::Info, true > LogVerbatim
HFWedgeSD(const std::string &, const SensitiveDetectorCatalog &clg, const SimTrackManager *)
Definition: HFWedgeSD.cc:23
std::vector< PCaloHit > PCaloHitContainer
int previousID
Definition: HFWedgeSD.h:46
void updateHit(HFShowerG4Hit *)
Definition: HFWedgeSD.cc:109
void Initialize(G4HCofThisEvent *HCE) override
Definition: HFWedgeSD.cc:30
int hcID
Definition: HFWedgeSD.h:42
G4ThreeVector momDir
Definition: HFWedgeSD.h:48
bool ProcessHits(G4Step *step, G4TouchableHistory *tHistory) override
Definition: HFWedgeSD.cc:42
void setLocalPos(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:42
void DrawAll() override
Definition: HFWedgeSD.cc:68
void fillHits(edm::PCaloHitContainer &, const std::string &) override
Definition: HFWedgeSD.cc:130
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:120
G4bool hitExists()
Definition: HFWedgeSD.cc:72
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:61
void setPrimMomDir(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:44
uint32_t setDetUnitId(const G4Step *) override
Definition: HFWedgeSD.cc:125
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:70
~HFWedgeSD() override
Definition: HFWedgeSD.cc:28
HFShowerG4Hit * createNewHit()
Definition: HFWedgeSD.cc:88
void clear() override
Definition: HFWedgeSD.cc:66
int currentID
Definition: HFWedgeSD.h:46
G4THitsCollection< HFShowerG4Hit > HFShowerG4HitsCollection
Definition: HFShowerG4Hit.h:55
G4ThreeVector globalPos
Definition: HFWedgeSD.h:48