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  const SensitiveDetectorCatalog& clg,
23  const SimTrackManager* manager)
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 << ":" << theHC;
37 
38  clearHits();
39 }
40 
41 G4bool HFWedgeSD::ProcessHits(G4Step* aStep, G4TouchableHistory*) {
42  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
43  const G4VTouchable* touch = preStepPoint->GetTouchable();
44  currentID = setDetUnitId(aStep);
45  trackID = aStep->GetTrack()->GetTrackID();
46  edep = aStep->GetTotalEnergyDeposit();
47  time = (preStepPoint->GetGlobalTime()) / ns;
48 
49  globalPos = preStepPoint->GetPosition();
50  localPos = touch->GetHistory()->GetTopTransform().TransformPoint(globalPos);
51  const G4DynamicParticle* particle = aStep->GetTrack()->GetDynamicParticle();
52  momDir = particle->GetMomentumDirection();
53 
54  if (hitExists() == false && edep > 0.)
56 
57  return true;
58 }
59 
60 void HFWedgeSD::EndOfEvent(G4HCofThisEvent* HCE) {
61  edm::LogVerbatim("FiberSim") << "HFWedgeSD: Sees" << theHC->entries() << " hits";
62  clear();
63 }
64 
66 
68 
70 
72  // Update if in the same detector, time-slice and for same track
73  if (currentID == previousID) {
75  return true;
76  }
77 
78  std::map<int, HFShowerG4Hit*>::const_iterator it = hitMap.find(currentID);
79  if (it != hitMap.end()) {
81  return true;
82  }
83 
84  return false;
85 }
86 
88 #ifdef EDM_ML_DEBUG
89  edm::LogVerbatim("FiberSim") << "HFWedgeSD::CreateNewHit for ID " << currentID << " Track " << trackID
90  << " Edep: " << edep / CLHEP::MeV << " MeV; Time: " << time << " ns; Position (local) "
91  << localPos << " (global ) " << globalPos << " direction " << momDir;
92 #endif
93  HFShowerG4Hit* aHit = new HFShowerG4Hit;
94  aHit->setHitId(currentID);
95  aHit->setTrackId(trackID);
96  aHit->setTime(time);
97  aHit->setLocalPos(localPos);
98  aHit->setGlobalPos(globalPos);
99  aHit->setPrimMomDir(momDir);
100  updateHit(aHit);
101 
102  theHC->insert(aHit);
103  hitMap.insert(std::pair<int, HFShowerG4Hit*>(previousID, aHit));
104 
105  return aHit;
106 }
107 
109  if (edep != 0) {
110  aHit->updateEnergy(edep);
111 #ifdef EDM_ML_DEBUG
112  edm::LogVerbatim("FiberSim") << "HFWedgeSD: Add energy deposit in " << currentID << " edep " << edep / CLHEP::MeV << " 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 
HFWedgeSD::hitMap
std::map< int, HFShowerG4Hit * > hitMap
Definition: HFWedgeSD.h:46
HFShowerG4Hit::setPrimMomDir
void setPrimMomDir(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:44
SimTrackManager
Definition: SimTrackManager.h:35
HFWedgeSD::edep
double edep
Definition: HFWedgeSD.h:49
HFWedgeSD::DrawAll
void DrawAll() override
Definition: HFWedgeSD.cc:67
HFWedgeSD::createNewHit
HFShowerG4Hit * createNewHit()
Definition: HFWedgeSD.cc:87
HFShowerG4HitsCollection
G4THitsCollection< HFShowerG4Hit > HFShowerG4HitsCollection
Definition: HFShowerG4Hit.h:55
HFWedgeSD.h
HFWedgeSD::hcID
int hcID
Definition: HFWedgeSD.h:44
HFShowerG4Hit::setLocalPos
void setLocalPos(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:42
HFWedgeSD::hitExists
G4bool hitExists()
Definition: HFWedgeSD.cc:71
MeV
const double MeV
HFShowerG4Hit
Definition: HFShowerG4Hit.h:15
SensitiveCaloDetector
Definition: SensitiveCaloDetector.h:10
HFWedgeSD::PrintAll
void PrintAll() override
Definition: HFWedgeSD.cc:69
HFWedgeSD::setDetUnitId
uint32_t setDetUnitId(const G4Step *) override
Definition: HFWedgeSD.cc:123
SensitiveDetectorCatalog
Definition: SensitiveDetectorCatalog.h:10
HFWedgeSD::HFWedgeSD
HFWedgeSD(const std::string &, const SensitiveDetectorCatalog &clg, const SimTrackManager *)
Definition: HFWedgeSD.cc:21
HFWedgeSD::momDir
G4ThreeVector momDir
Definition: HFWedgeSD.h:50
bysipixelclustmulteventfilter_cfi.collectionName
collectionName
Definition: bysipixelclustmulteventfilter_cfi.py:5
HFWedgeSD::time
double time
Definition: HFWedgeSD.h:49
HFWedgeSD::fillHits
void fillHits(edm::PCaloHitContainer &, const std::string &) override
Definition: HFWedgeSD.cc:128
HFWedgeSD::~HFWedgeSD
~HFWedgeSD() override
Definition: HFWedgeSD.cc:28
HFWedgeSD::currentHit
HFShowerG4Hit * currentHit
Definition: HFWedgeSD.h:51
HFWedgeSD::updateHit
void updateHit(HFShowerG4Hit *)
Definition: HFWedgeSD.cc:108
HFShowerG4Hit::setGlobalPos
void setGlobalPos(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:43
HFShowerG4Hit::setTrackId
void setTrackId(G4int trackId)
Definition: HFShowerG4Hit.h:38
HFShowerG4Hit::setTime
void setTime(G4double t)
Definition: HFShowerG4Hit.h:41
HFWedgeSD::clear
void clear() override
Definition: HFWedgeSD.cc:65
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
HFWedgeSD::clearHits
void clearHits() override
Definition: HFWedgeSD.cc:118
HFWedgeSD::globalPos
G4ThreeVector globalPos
Definition: HFWedgeSD.h:50
HFWedgeSD::ProcessHits
bool ProcessHits(G4Step *step, G4TouchableHistory *tHistory) override
Definition: HFWedgeSD.cc:41
HFWedgeSD::localPos
G4ThreeVector localPos
Definition: HFWedgeSD.h:50
HFShowerG4Hit::setHitId
void setHitId(G4int hitId)
Definition: HFShowerG4Hit.h:37
edm::LogVerbatim
Log< level::Info, true > LogVerbatim
Definition: MessageLogger.h:128
HFWedgeSD::Initialize
void Initialize(G4HCofThisEvent *HCE) override
Definition: HFWedgeSD.cc:30
HFWedgeSD::currentID
int currentID
Definition: HFWedgeSD.h:48
Point3D.h
edm::PCaloHitContainer
std::vector< PCaloHit > PCaloHitContainer
Definition: PCaloHitContainer.h:8
HFWedgeSD::EndOfEvent
void EndOfEvent(G4HCofThisEvent *eventHC) override
Definition: HFWedgeSD.cc:60
HFWedgeSD::trackID
int trackID
Definition: HFWedgeSD.h:48
HFWedgeSD::previousID
int previousID
Definition: HFWedgeSD.h:48
HFWedgeSD::theHC
HFShowerG4HitsCollection * theHC
Definition: HFWedgeSD.h:45
HFShowerG4Hit::updateEnergy
void updateEnergy(G4double edep)
Definition: HFShowerG4Hit.h:40