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 
20  const edm::EventSetup& es,
21  const SensitiveDetectorCatalog& clg,
22  edm::ParameterSet const& p,
23  const SimTrackManager* manager)
24  : SensitiveCaloDetector(iname, es, clg, p),
25  m_trackManager(manager),
26  hcID(-1),
27  theHC(nullptr),
28  currentHit(nullptr) {}
29 
31 
32 void HFWedgeSD::Initialize(G4HCofThisEvent* HCE) {
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  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  LogDebug("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  LogDebug("FiberSim") << "HFWedgeSD::CreateNewHit for ID " << currentID << " Track " << trackID
90  << " Edep: " << edep / MeV << " MeV; Time: " << time << " ns; Position (local) " << localPos
91  << " (global ) " << globalPos << " direction " << momDir;
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  LogDebug("FiberSim") << "HFWedgeSD: Add energy deposit in " << currentID << " edep " << edep / MeV << " MeV";
111  }
113 }
114 
116  hitMap.erase(hitMap.begin(), hitMap.end());
117  previousID = -1;
118 }
119 
120 uint32_t HFWedgeSD::setDetUnitId(const G4Step* aStep) {
121  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
122  return (touch->GetReplicaNumber(0));
123 }
124 
HFWedgeSD::hitMap
std::map< int, HFShowerG4Hit * > hitMap
Definition: HFWedgeSD.h:48
HFShowerG4Hit::setPrimMomDir
void setPrimMomDir(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:44
SimTrackManager
Definition: SimTrackManager.h:35
HFWedgeSD::edep
double edep
Definition: HFWedgeSD.h:51
HFWedgeSD::DrawAll
void DrawAll() override
Definition: HFWedgeSD.cc:68
HFWedgeSD::createNewHit
HFShowerG4Hit * createNewHit()
Definition: HFWedgeSD.cc:88
HFShowerG4HitsCollection
G4THitsCollection< HFShowerG4Hit > HFShowerG4HitsCollection
Definition: HFShowerG4Hit.h:55
HFWedgeSD.h
HFWedgeSD::hcID
int hcID
Definition: HFWedgeSD.h:46
AlCaHLTBitMon_ParallelJobs.p
p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
HFShowerG4Hit::setLocalPos
void setLocalPos(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:42
HFWedgeSD::hitExists
G4bool hitExists()
Definition: HFWedgeSD.cc:72
MeV
const double MeV
HFShowerG4Hit
Definition: HFShowerG4Hit.h:15
SensitiveCaloDetector
Definition: SensitiveCaloDetector.h:10
HFWedgeSD::PrintAll
void PrintAll() override
Definition: HFWedgeSD.cc:70
HFWedgeSD::setDetUnitId
uint32_t setDetUnitId(const G4Step *) override
Definition: HFWedgeSD.cc:120
SensitiveDetectorCatalog
Definition: SensitiveDetectorCatalog.h:10
HFWedgeSD::momDir
G4ThreeVector momDir
Definition: HFWedgeSD.h:52
bysipixelclustmulteventfilter_cfi.collectionName
collectionName
Definition: bysipixelclustmulteventfilter_cfi.py:5
HFWedgeSD::time
double time
Definition: HFWedgeSD.h:51
HFWedgeSD::fillHits
void fillHits(edm::PCaloHitContainer &, const std::string &) override
Definition: HFWedgeSD.cc:125
AlCaHLTBitMon_QueryRunRegistry.string
string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
HFWedgeSD::~HFWedgeSD
~HFWedgeSD() override
Definition: HFWedgeSD.cc:30
HFWedgeSD::currentHit
HFShowerG4Hit * currentHit
Definition: HFWedgeSD.h:53
HFWedgeSD::updateHit
void updateHit(HFShowerG4Hit *)
Definition: HFWedgeSD.cc:107
LogDebug
#define LogDebug(id)
Definition: MessageLogger.h:670
HFShowerG4Hit::setGlobalPos
void setGlobalPos(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:43
edm::ParameterSet
Definition: ParameterSet.h:36
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:66
edm::EventSetup
Definition: EventSetup.h:57
HFWedgeSD::clearHits
void clearHits() override
Definition: HFWedgeSD.cc:115
HFWedgeSD::globalPos
G4ThreeVector globalPos
Definition: HFWedgeSD.h:52
HFWedgeSD::ProcessHits
bool ProcessHits(G4Step *step, G4TouchableHistory *tHistory) override
Definition: HFWedgeSD.cc:42
HFWedgeSD::localPos
G4ThreeVector localPos
Definition: HFWedgeSD.h:52
HFShowerG4Hit::setHitId
void setHitId(G4int hitId)
Definition: HFShowerG4Hit.h:37
HFWedgeSD::Initialize
void Initialize(G4HCofThisEvent *HCE) override
Definition: HFWedgeSD.cc:32
HFWedgeSD::currentID
int currentID
Definition: HFWedgeSD.h:50
Point3D.h
edm::PCaloHitContainer
std::vector< PCaloHit > PCaloHitContainer
Definition: PCaloHitContainer.h:8
HFWedgeSD::EndOfEvent
void EndOfEvent(G4HCofThisEvent *eventHC) override
Definition: HFWedgeSD.cc:61
HFWedgeSD::trackID
int trackID
Definition: HFWedgeSD.h:50
HFWedgeSD::HFWedgeSD
HFWedgeSD(const std::string &, const edm::EventSetup &cpv, const SensitiveDetectorCatalog &clg, edm::ParameterSet const &p, const SimTrackManager *)
Definition: HFWedgeSD.cc:19
HFWedgeSD::previousID
int previousID
Definition: HFWedgeSD.h:50
HFWedgeSD::theHC
HFShowerG4HitsCollection * theHC
Definition: HFWedgeSD.h:47
HFShowerG4Hit::updateEnergy
void updateEnergy(G4double edep)
Definition: HFShowerG4Hit.h:40