CMS 3D CMS Logo

PPSPixelSD.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: PPS
4 // Class : PPSPixelSD
5 //
6 // Implementation:
7 // <Notes on implementation>
8 //
9 // Original Author: F.Ferro
10 // Created: May 4, 2015
11 //
12 
13 // system include files
14 
15 // user include files
18 
22 
25 
28 
29 #include "G4SDManager.hh"
30 #include "G4Step.hh"
31 #include "G4Track.hh"
32 #include "G4VProcess.hh"
33 
34 #include "G4PhysicalConstants.hh"
35 #include "G4SystemOfUnits.hh"
36 
38  const SensitiveDetectorCatalog& clg,
39  edm::ParameterSet const& p,
40  SimTrackManager const* manager)
41  : SensitiveTkDetector(pname, clg), theManager_(manager) {
42  //Add PPS Sentitive Detector Names
43  collectionName.insert(pname);
44 
45  //Parameters
46  edm::ParameterSet m_p = p.getParameter<edm::ParameterSet>("PPSPixelSD");
47  int verbn = m_p.getUntrackedParameter<int>("Verbosity");
48  SetVerboseLevel(verbn);
49  slave_ = std::make_unique<TrackingSlaveSD>(pname);
50  if (pname == "CTPPSPixelHits") {
51  numberingScheme_ = std::make_unique<PPSPixelOrganization>();
52  } else {
53  edm::LogWarning("PPSSim") << "PPSPixelSD: ReadoutName " << pname << " not supported";
54  }
55  edm::LogVerbatim("PPSSim") << "PPSPixelSD: Instantiation completed for " << pname;
56 }
57 
59 
60 bool PPSPixelSD::ProcessHits(G4Step* aStep, G4TouchableHistory*) {
61  eloss_ = aStep->GetTotalEnergyDeposit();
62  if (eloss_ > 0.0) {
63  eloss_ /= GeV;
64  stepInfo(aStep);
65  LogDebug("PPSSim") << "PPSPixelSD: ProcessHits: edep=" << eloss_ << " "
66  << theTrack_->GetDefinition()->GetParticleName();
67  if (!hitExists())
68  createNewHit();
69  }
70  return true;
71 }
72 
73 uint32_t PPSPixelSD::setDetUnitId(const G4Step* aStep) {
74  return (numberingScheme_ == nullptr ? 0 : numberingScheme_->unitID(aStep));
75 }
76 
77 void PPSPixelSD::Initialize(G4HCofThisEvent* HCE) {
78  LogDebug("PPSSim") << "PPSPixelSD : Initialize called for " << GetName();
79 
80  theHC_ = new PPSPixelG4HitCollection(GetName(), collectionName[0]);
81  G4SDManager::GetSDMpointer()->AddNewCollection(GetName(), collectionName[0]);
82 
83  if (hcID_ < 0)
84  hcID_ = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
85  HCE->AddHitsCollection(hcID_, theHC_);
86 
87  tsID_ = -2;
88  edm::LogVerbatim("PPSSim") << "PPSPixelSD: is initialized " << GetName();
89 }
90 
91 void PPSPixelSD::EndOfEvent(G4HCofThisEvent*) {
92  // here we loop over transient hits and make them persistent
93  for (unsigned int j = 0; j < (unsigned int)theHC_->entries(); ++j) {
94  PPSPixelG4Hit* aHit = (*theHC_)[j];
95 #ifdef debug
96  LogDebug("PPSSim") << "HIT NUMERO " << j << "unit ID = " << aHit->unitID() << "\n"
97  << " "
98  << "enrty z " << aHit->entry().z() << "\n"
99  << " "
100  << "theta " << aHit->thetaAtEntry() << "\n";
101 #endif
102 
103  Local3DPoint Enter(aHit->entryPoint().x(), aHit->entryPoint().y(), aHit->entryPoint().z());
104  Local3DPoint Exit(aHit->exitPoint().x(), aHit->exitPoint().y(), aHit->exitPoint().z());
105  slave_->processHits(PSimHit(Enter,
106  Exit,
107  aHit->p(),
108  aHit->tof(),
109  aHit->energyLoss(),
110  aHit->particleType(),
111  aHit->unitID(),
112  aHit->trackID(),
113  aHit->thetaAtEntry(),
114  aHit->phiAtEntry()));
115  }
116 }
117 
119  LogDebug("PPSSim") << "PPSPixelSD: Collection " << theHC_->GetName();
120  theHC_->PrintAllHits();
121 }
122 
124  if (slave_->name() == n) {
125  c = slave_->hits();
126  }
127 }
128 
130  LogDebug("PPSSim") << "PPSPixelSD: dispatched BeginOfEvent for " << GetName();
131  clearHits();
132  eventno_ = (*i)()->GetEventID();
133 }
134 
135 void PPSPixelSD::update(const ::EndOfEvent*) {}
136 
137 void PPSPixelSD::clearHits() { slave_->Initialize(); }
138 
139 G4ThreeVector PPSPixelSD::setToLocal(const G4ThreeVector& global) {
140  G4ThreeVector localPoint;
141  const G4VTouchable* touch = preStepPoint_->GetTouchable();
142  localPoint = touch->GetHistory()->GetTopTransform().TransformPoint(global);
143  return localPoint;
144 }
145 
146 void PPSPixelSD::stepInfo(const G4Step* aStep) {
147  preStepPoint_ = aStep->GetPreStepPoint();
148  postStepPoint_ = aStep->GetPostStepPoint();
149  theTrack_ = aStep->GetTrack();
150 
151  hitPoint_ = preStepPoint_->GetPosition();
152  exitPoint_ = postStepPoint_->GetPosition();
153  currentPV_ = preStepPoint_->GetPhysicalVolume();
156 
157  tSlice_ = postStepPoint_->GetGlobalTime() / nanosecond;
158  tSliceID_ = (int)tSlice_;
159  unitID_ = setDetUnitId(aStep);
160 #ifdef debug
161  LogDebug("PPSSim") << "UNIT ID " << unitID_;
162 #endif
163  primaryID_ = theTrack_->GetTrackID();
164  parentID_ = theTrack_->GetParentID();
165 
166  incidentEnergy_ = theTrack_->GetTotalEnergy() / GeV;
167 
168  pabs_ = aStep->GetPreStepPoint()->GetMomentum().mag() / GeV;
169  tof_ = aStep->GetPostStepPoint()->GetGlobalTime() / nanosecond;
170 
171  particleType_ = theTrack_->GetDefinition()->GetPDGEncoding();
172 
173  thetaAtEntry_ = aStep->GetPreStepPoint()->GetPosition().theta();
174  phiAtEntry_ = aStep->GetPreStepPoint()->GetPosition().phi();
175 
176  vx_ = theTrack_->GetVertexPosition().x();
177  vy_ = theTrack_->GetVertexPosition().y();
178  vz_ = theTrack_->GetVertexPosition().z();
179 }
180 
182  // Update if in the same detector, time-slice and for same track
183  if (tSliceID_ == tsID_ && unitID_ == previousUnitID_) {
184  updateHit();
185  return true;
186  }
187 
188  // look in the HitContainer whether a hit with the same unitID_,
189  // tSliceID_ already exists: secondary energy deposition
190  // is added to the primary hit
191  int nhits = theHC_->entries();
192  for (int j = 0; j < nhits; ++j) {
193  PPSPixelG4Hit* aPreviousHit = (*theHC_)[j];
194  if (aPreviousHit->timeSliceID() == tSliceID_ && aPreviousHit->unitID() == unitID_) {
195  currentHit_ = aPreviousHit;
196  updateHit();
197  return true;
198  }
199  }
200  return false;
201 }
202 
204 #ifdef debug
205  LogDebug("PPSSim") << "PPSPixelSD::CreateNewHit for"
206  << " PV " << currentPV_->GetName() << " PVid = " << currentPV_->GetCopyNo()
207  << " MVid = " << currentPV_->GetMother()->GetCopyNo() << " Unit " << unitID_ << "\n"
208  << " primary " << primaryID_ << " time slice " << tSliceID_ << " For Track "
209  << " which is a " << theTrack_->GetDefinition()->GetParticleName();
210 
211  if (parentID == 0) {
212  LogDebug("PPSSim") << "PPSPixelSD: primary of energy " << incidentEnergy_;
213  } else {
214  LogDebug("PPSSim") << " daughter of track: " << parentID;
215  }
216 
217  if (theTrack_->GetCreatorProcess() != nullptr)
218  LogDebug("PPSSim") << theTrack_->GetCreatorProcess()->GetProcessName();
219  else
220  LogDebug("PPSSim") << "NO process";
221 #endif
222 
228 
234 
235  G4ThreeVector pos = 0.5 * (hitPoint_ + exitPoint_);
239 
244 
245  updateHit();
247 }
248 
250 #ifdef debug
251  LogDebug("PPSSim") << "G4PPSPixelSD updateHit: add eloss=" << eloss_ << "\nCurrentHit=" << currentHit_
252  << ", PostStepPoint=" << postStepPoint_->GetPosition();
253 #endif
255  // buffer for next steps:
256  tsID_ = tSliceID_;
258 }
259 
261  if (hit == nullptr) {
262  edm::LogWarning("PPSSim") << "PPSPixelSD: hit to be stored is NULL !!";
263  return;
264  }
265 
266  theHC_->insert(hit);
267 }
void setParticleType(short i)
float z() const
Log< level::Info, true > LogVerbatim
PPSPixelSD(const std::string &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, SimTrackManager const *)
Definition: PPSPixelSD.cc:37
uint32_t previousUnitID_
Definition: PPSPixelSD.h:107
G4ThreeVector setToLocal(const G4ThreeVector &globalPoint)
Definition: PPSPixelSD.cc:139
uint32_t setDetUnitId(const G4Step *) override
Definition: PPSPixelSD.cc:73
void stepInfo(const G4Step *aStep)
Definition: PPSPixelSD.cc:146
int parentID_
Definition: PPSPixelSD.h:110
G4ThreeVector exitPoint_
Definition: PPSPixelSD.h:90
G4int hcID_
Definition: PPSPixelSD.h:85
void clearHits() override
Definition: PPSPixelSD.cc:137
void PrintAll() override
Definition: PPSPixelSD.cc:118
void setExitPoint(const G4ThreeVector &)
float energyLoss() const
G4VPhysicalVolume * currentPV_
Definition: PPSPixelSD.h:84
G4StepPoint * preStepPoint_
Definition: PPSPixelSD.h:87
uint32_t unitID() const
float vz_
Definition: PPSPixelSD.h:104
int particleType() const
float phiAtEntry_
Definition: PPSPixelSD.h:101
void setP(float e)
void createNewHit()
Definition: PPSPixelSD.cc:203
void storeHit(PPSPixelG4Hit *)
Definition: PPSPixelSD.cc:260
bool ProcessHits(G4Step *, G4TouchableHistory *) override
Definition: PPSPixelSD.cc:60
void setVy(float p)
T getUntrackedParameter(std::string const &, T const &) const
const G4ThreeVector & entryPoint() const
float tof_
Definition: PPSPixelSD.h:98
void update(const BeginOfEvent *) override
This routine will be called when the appropriate signal arrives.
Definition: PPSPixelSD.cc:129
float vx_
Definition: PPSPixelSD.h:102
const G4ThreeVector & exitPoint() const
int eventno_
Definition: PPSPixelSD.h:112
int primaryID_
Definition: PPSPixelSD.h:109
PPSPixelG4Hit * currentHit_
Definition: PPSPixelSD.h:82
void setPhiAtEntry(float f)
G4Track * theTrack_
Definition: PPSPixelSD.h:83
float vy_
Definition: PPSPixelSD.h:103
float thetaAtEntry_
Definition: PPSPixelSD.h:100
float p() const
G4ThreeVector theLocalExitPoint_
Definition: PPSPixelSD.h:92
float thetaAtEntry() const
void setParentId(int p)
bool hitExists()
Definition: PPSPixelSD.cc:181
void Initialize(G4HCofThisEvent *HCE) override
Definition: PPSPixelSD.cc:77
float incidentEnergy_
Definition: PPSPixelSD.h:96
std::unique_ptr< TrackingSlaveSD > slave_
Definition: PPSPixelSD.h:77
void setTof(float e)
int timeSliceID() const
float pabs_
Definition: PPSPixelSD.h:97
void setTimeSlice(double d)
void setEntryPoint(const G4ThreeVector &)
PPSPixelG4HitCollection * theHC_
Definition: PPSPixelSD.h:80
double tSlice_
Definition: PPSPixelSD.h:94
float tof() const
void setThetaAtEntry(float t)
std::unique_ptr< PPSVDetectorOrganization > numberingScheme_
Definition: PPSPixelSD.h:78
void setVx(float p)
void fillHits(edm::PSimHitContainer &, const std::string &) override
Definition: PPSPixelSD.cc:123
~PPSPixelSD() override
Definition: PPSPixelSD.cc:58
void setTrackID(int i)
int tSliceID_
Definition: PPSPixelSD.h:111
void setUnitID(uint32_t i)
G4ThreeVector hitPoint_
Definition: PPSPixelSD.h:89
G4THitsCollection< PPSPixelG4Hit > PPSPixelG4HitCollection
void setVz(float p)
G4StepPoint * postStepPoint_
Definition: PPSPixelSD.h:88
short particleType_
Definition: PPSPixelSD.h:113
std::vector< PSimHit > PSimHitContainer
uint32_t unitID_
Definition: PPSPixelSD.h:106
void EndOfEvent(G4HCofThisEvent *eventHC) override
Definition: PPSPixelSD.cc:91
int trackID() const
Log< level::Warning, false > LogWarning
void updateHit()
Definition: PPSPixelSD.cc:249
void setIncidentEnergy(double e)
void setEnergyLoss(float e)
G4ThreeVector theLocalEntryPoint_
Definition: PPSPixelSD.h:91
void setMeanPosition(const G4ThreeVector &a)
Definition: PPSPixelG4Hit.h:54
double eloss_
Definition: PPSPixelSD.h:95
#define LogDebug(id)
float phiAtEntry() const