CMS 3D CMS Logo

HFChamberSD.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(name, clg), theHCID(-1), theHC(nullptr), theNSteps(0) {
25  edm::LogVerbatim("FiberSim") << "HFChamberSD : Instantiated for " << name;
26 }
27 
29 
30 void HFChamberSD::Initialize(G4HCofThisEvent* HCE) {
31  edm::LogVerbatim("FiberSim") << "HFChamberSD : Initialize called for " << GetName() << " in collection " << HCE;
32  theHC = new HFShowerG4HitsCollection(GetName(), collectionName[0]);
33  if (theHCID < 0)
34  theHCID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
35  HCE->AddHitsCollection(theHCID, theHC);
36  edm::LogVerbatim("FiberSim") << "HFChamberSD : Add hit collectrion for " << collectionName[0] << ":" << theHCID << ":" << theHC;
37 }
38 
39 G4bool HFChamberSD::ProcessHits(G4Step* aStep, G4TouchableHistory*) {
40  //do not process hits other than primary particle hits:
41  double charge = aStep->GetTrack()->GetDefinition()->GetPDGCharge();
42  int trackID = aStep->GetTrack()->GetTrackID();
43  if (charge == 0. || trackID != 1 || aStep->GetTrack()->GetParentID() != 0 ||
44  aStep->GetTrack()->GetCreatorProcess() != nullptr)
45  return false;
46  ++theNSteps;
47  //if(theNSteps>1)return false;
48 
49  G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
50  const G4VTouchable* touch = preStepPoint->GetTouchable();
51  int detID = setDetUnitId(aStep);
52 
53  double edep = aStep->GetTotalEnergyDeposit();
54  double time = (preStepPoint->GetGlobalTime()) / ns;
55 
56  const G4ThreeVector& globalPos = preStepPoint->GetPosition();
57  G4ThreeVector localPos = touch->GetHistory()->GetTopTransform().TransformPoint(globalPos);
58  const G4DynamicParticle* particle = aStep->GetTrack()->GetDynamicParticle();
59  const G4ThreeVector& momDir = particle->GetMomentumDirection();
60 
61  HFShowerG4Hit* aHit = new HFShowerG4Hit(detID, trackID, edep, time);
62  aHit->setLocalPos(localPos);
63  aHit->setGlobalPos(globalPos);
64  aHit->setPrimMomDir(momDir);
65 
66 #ifdef EDM_ML_DEBUG
67  edm::LogVerbatim("FiberSim") << "HFChamberSD: Hit created in (" << touch->GetVolume(0)->GetLogicalVolume()->GetName() << ") ID " << detID << " Track " << trackID << " Edep: " << edep / CLHEP::MeV << " MeV; Time: " << time << " ns; Position (local) " << localPos << " (global ) " << globalPos << " direction " << momDir;
68 #endif
69  theHC->insert(aHit);
70  return true;
71 }
72 
73 void HFChamberSD::EndOfEvent(G4HCofThisEvent* HCE) {
74  edm::LogVerbatim("FiberSim") << "HFChamberSD: Finds " << theHC->entries() << " hits";
75  clear();
76 }
77 
79 
81 
83 
85 
86 uint32_t HFChamberSD::setDetUnitId(const G4Step* aStep) {
87  const G4VTouchable* touch = aStep->GetPreStepPoint()->GetTouchable();
88  return (touch->GetReplicaNumber(0));
89 }
90 
HFChamberSD::clearHits
void clearHits() override
Definition: HFChamberSD.cc:84
HFShowerG4Hit::setPrimMomDir
void setPrimMomDir(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:44
SimTrackManager
Definition: SimTrackManager.h:35
HFShowerG4HitsCollection
G4THitsCollection< HFShowerG4Hit > HFShowerG4HitsCollection
Definition: HFShowerG4Hit.h:55
HFShowerG4Hit::setLocalPos
void setLocalPos(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:42
protons_cff.time
time
Definition: protons_cff.py:35
MeV
const double MeV
HFChamberSD::HFChamberSD
HFChamberSD(const std::string &, const SensitiveDetectorCatalog &, const SimTrackManager *)
Definition: HFChamberSD.cc:21
HFShowerG4Hit
Definition: HFShowerG4Hit.h:15
SensitiveCaloDetector
Definition: SensitiveCaloDetector.h:10
HFChamberSD::~HFChamberSD
~HFChamberSD() override
Definition: HFChamberSD.cc:28
HFChamberSD::Initialize
void Initialize(G4HCofThisEvent *HCE) override
Definition: HFChamberSD.cc:30
SensitiveDetectorCatalog
Definition: SensitiveDetectorCatalog.h:10
HFChamberSD::theHC
HFShowerG4HitsCollection * theHC
Definition: HFChamberSD.h:41
bysipixelclustmulteventfilter_cfi.collectionName
collectionName
Definition: bysipixelclustmulteventfilter_cfi.py:5
ALCARECOTkAlJpsiMuMu_cff.charge
charge
Definition: ALCARECOTkAlJpsiMuMu_cff.py:47
HFChamberSD::setDetUnitId
uint32_t setDetUnitId(const G4Step *) override
Definition: HFChamberSD.cc:86
HFChamberSD::fillHits
void fillHits(edm::PCaloHitContainer &, const std::string &) override
Definition: HFChamberSD.cc:91
HFChamberSD::theHCID
G4int theHCID
Definition: HFChamberSD.h:40
HFShowerG4Hit::setGlobalPos
void setGlobalPos(const G4ThreeVector &xyz)
Definition: HFShowerG4Hit.h:43
HFChamberSD.h
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
HFChamberSD::theNSteps
int theNSteps
Definition: HFChamberSD.h:42
HFChamberSD::DrawAll
void DrawAll() override
Definition: HFChamberSD.cc:80
edm::LogVerbatim
Log< level::Info, true > LogVerbatim
Definition: MessageLogger.h:128
HFChamberSD::PrintAll
void PrintAll() override
Definition: HFChamberSD.cc:82
Point3D.h
edm::PCaloHitContainer
std::vector< PCaloHit > PCaloHitContainer
Definition: PCaloHitContainer.h:8
Skims_PA_cff.name
name
Definition: Skims_PA_cff.py:17
HFChamberSD::ProcessHits
G4bool ProcessHits(G4Step *aStep, G4TouchableHistory *ROhist) override
Definition: HFChamberSD.cc:39
HFChamberSD::EndOfEvent
void EndOfEvent(G4HCofThisEvent *HCE) override
Definition: HFChamberSD.cc:73
HFChamberSD::clear
void clear() override
Definition: HFChamberSD.cc:78