CMS 3D CMS Logo

List of all members | Public Member Functions | Protected Member Functions | Private Member Functions | Private Attributes
PPSPixelSD Class Reference

#include <SimG4CMS/PPS/interface/PPSPixelSD.h>

Inheritance diagram for PPSPixelSD:
SensitiveTkDetector Observer< const BeginOfEvent *> Observer< const EndOfEvent *> SensitiveDetector

Public Member Functions

void clearHits () override
 
void EndOfEvent (G4HCofThisEvent *eventHC) override
 
void fillHits (edm::PSimHitContainer &, const std::string &) override
 
void Initialize (G4HCofThisEvent *HCE) override
 
 PPSPixelSD (const std::string &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, SimTrackManager const *)
 
void PrintAll () override
 
bool ProcessHits (G4Step *, G4TouchableHistory *) override
 
uint32_t setDetUnitId (const G4Step *) override
 
 ~PPSPixelSD () override
 
- Public Member Functions inherited from SensitiveTkDetector
 SensitiveTkDetector (const std::string &iname, const SensitiveDetectorCatalog &clg)
 
- Public Member Functions inherited from SensitiveDetector
void EndOfEvent (G4HCofThisEvent *eventHC) override
 
const std::vector< std::string > & getNames () const
 
void Initialize (G4HCofThisEvent *eventHC) override
 
bool isCaloSD () const
 
 SensitiveDetector (const std::string &iname, const SensitiveDetectorCatalog &, bool calo, const std::string &newcollname="")
 
 ~SensitiveDetector () override
 
- Public Member Functions inherited from Observer< const BeginOfEvent *>
 Observer ()
 
void slotForUpdate (const BeginOfEvent * iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const EndOfEvent *>
 Observer ()
 
void slotForUpdate (const EndOfEvent * iT)
 
virtual ~Observer ()
 

Protected Member Functions

void update (const BeginOfEvent *) override
 This routine will be called when the appropriate signal arrives. More...
 
void update (const ::EndOfEvent *) override
 
- Protected Member Functions inherited from SensitiveDetector
TrackInformationcmsTrackInformation (const G4Track *aTrack)
 
Local3DPoint ConvertToLocal3DPoint (const G4ThreeVector &point) const
 
Local3DPoint FinalStepPosition (const G4Step *step, coordinates) const
 
Local3DPoint InitialStepPosition (const G4Step *step, coordinates) const
 
Local3DPoint LocalPostStepPosition (const G4Step *step) const
 
Local3DPoint LocalPreStepPosition (const G4Step *step) const
 
void NaNTrap (const G4Step *step) const
 
void setNames (const std::vector< std::string > &)
 
- Protected Member Functions inherited from Observer< const EndOfEvent *>
virtual void update (const EndOfEvent *)=0
 This routine will be called when the appropriate signal arrives. More...
 

Private Member Functions

void createNewHit ()
 
bool hitExists ()
 
G4ThreeVector setToLocal (const G4ThreeVector &globalPoint)
 
void stepInfo (const G4Step *aStep)
 
void storeHit (PPSPixelG4Hit *)
 
void updateHit ()
 

Private Attributes

PPSPixelG4HitcurrentHit_ = nullptr
 
G4VPhysicalVolume * currentPV_ = nullptr
 
double eloss_ = 0.0
 
int eventno_ = 0
 
G4ThreeVector exitPoint_
 
G4int hcID_ = -1
 
G4ThreeVector hitPoint_
 
float incidentEnergy_ = 0.f
 
std::unique_ptr< PPSVDetectorOrganizationnumberingScheme_
 
float pabs_ = 0.f
 
int parentID_ = 0
 
short particleType_ = 0
 
float phiAtEntry_ = 0.f
 
G4StepPoint * postStepPoint_ = nullptr
 
G4StepPoint * preStepPoint_ = nullptr
 
uint32_t previousUnitID_ = 0
 
int primaryID_ = 0
 
std::unique_ptr< TrackingSlaveSDslave_
 
PPSPixelG4HitCollectiontheHC_ = nullptr
 
G4ThreeVector theLocalEntryPoint_
 
G4ThreeVector theLocalExitPoint_
 
const SimTrackManagertheManager_
 
float thetaAtEntry_ = 0.f
 
G4Track * theTrack_ = nullptr
 
float tof_ = 0.f
 
int tsID_ = 0
 
double tSlice_ = 0.0
 
int tSliceID_ = 0
 
uint32_t unitID_ = 0
 
float vx_ = 0.f
 
float vy_ = 0.f
 
float vz_ = 0.f
 

Additional Inherited Members

- Protected Types inherited from SensitiveDetector
enum  coordinates { WorldCoordinates, LocalCoordinates }
 

Detailed Description

Description: Stores hits of PPSPixel in appropriate container

Usage: Used in sensitive detector builder

Definition at line 46 of file PPSPixelSD.h.

Constructor & Destructor Documentation

◆ PPSPixelSD()

PPSPixelSD::PPSPixelSD ( const std::string &  pname,
const SensitiveDetectorCatalog clg,
edm::ParameterSet const &  p,
SimTrackManager const *  manager 
)

Definition at line 37 of file PPSPixelSD.cc.

References bysipixelclustmulteventfilter_cfi::collectionName, edm::ParameterSet::getUntrackedParameter(), numberingScheme_, AlCaHLTBitMon_ParallelJobs::p, unpackData-CaloStage2::pname, and slave_.

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 }
Log< level::Info, true > LogVerbatim
const SimTrackManager * theManager_
Definition: PPSPixelSD.h:81
T getUntrackedParameter(std::string const &, T const &) const
std::unique_ptr< TrackingSlaveSD > slave_
Definition: PPSPixelSD.h:77
std::unique_ptr< PPSVDetectorOrganization > numberingScheme_
Definition: PPSPixelSD.h:78
Log< level::Warning, false > LogWarning
SensitiveTkDetector(const std::string &iname, const SensitiveDetectorCatalog &clg)

◆ ~PPSPixelSD()

PPSPixelSD::~PPSPixelSD ( )
override

Definition at line 58 of file PPSPixelSD.cc.

58 {}

Member Function Documentation

◆ clearHits()

void PPSPixelSD::clearHits ( )
overridevirtual

Implements SensitiveDetector.

Definition at line 137 of file PPSPixelSD.cc.

References slave_.

Referenced by update().

137 { slave_->Initialize(); }
std::unique_ptr< TrackingSlaveSD > slave_
Definition: PPSPixelSD.h:77

◆ createNewHit()

void PPSPixelSD::createNewHit ( )
private

Definition at line 203 of file PPSPixelSD.cc.

References currentHit_, currentPV_, exitPoint_, hitPoint_, incidentEnergy_, LogDebug, pabs_, parentID_, particleType_, phiAtEntry_, primaryID_, PPSPixelG4Hit::setEntryPoint(), PPSPixelG4Hit::setExitPoint(), PPSPixelG4Hit::setIncidentEnergy(), PPSPixelG4Hit::setMeanPosition(), PPSPixelG4Hit::setP(), PPSPixelG4Hit::setParentId(), PPSPixelG4Hit::setParticleType(), PPSPixelG4Hit::setPhiAtEntry(), PPSPixelG4Hit::setThetaAtEntry(), PPSPixelG4Hit::setTimeSlice(), PPSPixelG4Hit::setTof(), PPSPixelG4Hit::setTrackID(), PPSPixelG4Hit::setUnitID(), PPSPixelG4Hit::setVx(), PPSPixelG4Hit::setVy(), PPSPixelG4Hit::setVz(), storeHit(), theLocalEntryPoint_, theLocalExitPoint_, thetaAtEntry_, theTrack_, tof_, tSlice_, tSliceID_, unitID_, updateHit(), vx_, vy_, and vz_.

Referenced by ProcessHits().

203  {
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 }
void setParticleType(short i)
int parentID_
Definition: PPSPixelSD.h:110
G4ThreeVector exitPoint_
Definition: PPSPixelSD.h:90
void setExitPoint(const G4ThreeVector &)
G4VPhysicalVolume * currentPV_
Definition: PPSPixelSD.h:84
float vz_
Definition: PPSPixelSD.h:104
float phiAtEntry_
Definition: PPSPixelSD.h:101
void setP(float e)
void storeHit(PPSPixelG4Hit *)
Definition: PPSPixelSD.cc:260
void setVy(float p)
float tof_
Definition: PPSPixelSD.h:98
float vx_
Definition: PPSPixelSD.h:102
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
G4ThreeVector theLocalExitPoint_
Definition: PPSPixelSD.h:92
void setParentId(int p)
float incidentEnergy_
Definition: PPSPixelSD.h:96
void setTof(float e)
float pabs_
Definition: PPSPixelSD.h:97
void setTimeSlice(double d)
void setEntryPoint(const G4ThreeVector &)
double tSlice_
Definition: PPSPixelSD.h:94
void setThetaAtEntry(float t)
void setVx(float p)
void setTrackID(int i)
int tSliceID_
Definition: PPSPixelSD.h:111
void setUnitID(uint32_t i)
G4ThreeVector hitPoint_
Definition: PPSPixelSD.h:89
void setVz(float p)
short particleType_
Definition: PPSPixelSD.h:113
uint32_t unitID_
Definition: PPSPixelSD.h:106
void updateHit()
Definition: PPSPixelSD.cc:249
void setIncidentEnergy(double e)
G4ThreeVector theLocalEntryPoint_
Definition: PPSPixelSD.h:91
void setMeanPosition(const G4ThreeVector &a)
Definition: PPSPixelG4Hit.h:54
#define LogDebug(id)

◆ EndOfEvent()

void PPSPixelSD::EndOfEvent ( G4HCofThisEvent *  eventHC)
override

Definition at line 91 of file PPSPixelSD.cc.

References PPSPixelG4Hit::energyLoss(), PPSPixelG4Hit::entryPoint(), PPSPixelG4Hit::exitPoint(), createfilelist::int, dqmiolumiharvest::j, LogDebug, PPSPixelG4Hit::p(), PPSPixelG4Hit::particleType(), PPSPixelG4Hit::phiAtEntry(), slave_, theHC_, PPSPixelG4Hit::thetaAtEntry(), PPSPixelG4Hit::tof(), PPSPixelG4Hit::trackID(), PPSPixelG4Hit::unitID(), and PPSPixelG4Hit::z().

91  {
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 }
float z() const
float energyLoss() const
uint32_t unitID() const
int particleType() const
const G4ThreeVector & entryPoint() const
const G4ThreeVector & exitPoint() const
float p() const
float thetaAtEntry() const
std::unique_ptr< TrackingSlaveSD > slave_
Definition: PPSPixelSD.h:77
PPSPixelG4HitCollection * theHC_
Definition: PPSPixelSD.h:80
float tof() const
int trackID() const
#define LogDebug(id)
float phiAtEntry() const

◆ fillHits()

void PPSPixelSD::fillHits ( edm::PSimHitContainer c,
const std::string &  n 
)
overridevirtual

Implements SensitiveTkDetector.

Definition at line 123 of file PPSPixelSD.cc.

References HltBtagPostValidation_cff::c, dqmiodumpmetadata::n, and slave_.

123  {
124  if (slave_->name() == n) {
125  c = slave_->hits();
126  }
127 }
std::unique_ptr< TrackingSlaveSD > slave_
Definition: PPSPixelSD.h:77

◆ hitExists()

bool PPSPixelSD::hitExists ( )
private

Definition at line 181 of file PPSPixelSD.cc.

References currentHit_, dqmiolumiharvest::j, nhits, previousUnitID_, theHC_, PPSPixelG4Hit::timeSliceID(), tsID_, tSliceID_, PPSPixelG4Hit::unitID(), unitID_, and updateHit().

Referenced by ProcessHits().

181  {
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 }
uint32_t previousUnitID_
Definition: PPSPixelSD.h:107
uint32_t unitID() const
PPSPixelG4Hit * currentHit_
Definition: PPSPixelSD.h:82
int timeSliceID() const
PPSPixelG4HitCollection * theHC_
Definition: PPSPixelSD.h:80
int tSliceID_
Definition: PPSPixelSD.h:111
uint32_t unitID_
Definition: PPSPixelSD.h:106
void updateHit()
Definition: PPSPixelSD.cc:249

◆ Initialize()

void PPSPixelSD::Initialize ( G4HCofThisEvent *  HCE)
override

Definition at line 77 of file PPSPixelSD.cc.

References bysipixelclustmulteventfilter_cfi::collectionName, hcID_, LogDebug, theHC_, and tsID_.

77  {
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 }
Log< level::Info, true > LogVerbatim
G4int hcID_
Definition: PPSPixelSD.h:85
PPSPixelG4HitCollection * theHC_
Definition: PPSPixelSD.h:80
G4THitsCollection< PPSPixelG4Hit > PPSPixelG4HitCollection
#define LogDebug(id)

◆ PrintAll()

void PPSPixelSD::PrintAll ( )
override

Definition at line 118 of file PPSPixelSD.cc.

References LogDebug, and theHC_.

118  {
119  LogDebug("PPSSim") << "PPSPixelSD: Collection " << theHC_->GetName();
120  theHC_->PrintAllHits();
121 }
PPSPixelG4HitCollection * theHC_
Definition: PPSPixelSD.h:80
#define LogDebug(id)

◆ ProcessHits()

bool PPSPixelSD::ProcessHits ( G4Step *  aStep,
G4TouchableHistory *   
)
overridevirtual

Implements SensitiveDetector.

Definition at line 60 of file PPSPixelSD.cc.

References createNewHit(), eloss_, hitExists(), LogDebug, stepInfo(), and theTrack_.

60  {
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 }
void stepInfo(const G4Step *aStep)
Definition: PPSPixelSD.cc:146
void createNewHit()
Definition: PPSPixelSD.cc:203
G4Track * theTrack_
Definition: PPSPixelSD.h:83
bool hitExists()
Definition: PPSPixelSD.cc:181
double eloss_
Definition: PPSPixelSD.h:95
#define LogDebug(id)

◆ setDetUnitId()

uint32_t PPSPixelSD::setDetUnitId ( const G4Step *  aStep)
overridevirtual

Implements SensitiveDetector.

Definition at line 73 of file PPSPixelSD.cc.

References numberingScheme_.

Referenced by stepInfo().

73  {
74  return (numberingScheme_ == nullptr ? 0 : numberingScheme_->unitID(aStep));
75 }
std::unique_ptr< PPSVDetectorOrganization > numberingScheme_
Definition: PPSPixelSD.h:78

◆ setToLocal()

G4ThreeVector PPSPixelSD::setToLocal ( const G4ThreeVector &  globalPoint)
private

Definition at line 139 of file PPSPixelSD.cc.

References preStepPoint_.

Referenced by stepInfo().

139  {
140  G4ThreeVector localPoint;
141  const G4VTouchable* touch = preStepPoint_->GetTouchable();
142  localPoint = touch->GetHistory()->GetTopTransform().TransformPoint(global);
143  return localPoint;
144 }
G4StepPoint * preStepPoint_
Definition: PPSPixelSD.h:87

◆ stepInfo()

void PPSPixelSD::stepInfo ( const G4Step *  aStep)
private

Definition at line 146 of file PPSPixelSD.cc.

References currentPV_, exitPoint_, hitPoint_, incidentEnergy_, createfilelist::int, LogDebug, pabs_, parentID_, particleType_, phiAtEntry_, postStepPoint_, preStepPoint_, primaryID_, setDetUnitId(), setToLocal(), theLocalEntryPoint_, theLocalExitPoint_, thetaAtEntry_, theTrack_, tof_, tSlice_, tSliceID_, unitID_, vx_, vy_, and vz_.

Referenced by ProcessHits().

146  {
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 }
G4ThreeVector setToLocal(const G4ThreeVector &globalPoint)
Definition: PPSPixelSD.cc:139
uint32_t setDetUnitId(const G4Step *) override
Definition: PPSPixelSD.cc:73
int parentID_
Definition: PPSPixelSD.h:110
G4ThreeVector exitPoint_
Definition: PPSPixelSD.h:90
G4VPhysicalVolume * currentPV_
Definition: PPSPixelSD.h:84
G4StepPoint * preStepPoint_
Definition: PPSPixelSD.h:87
float vz_
Definition: PPSPixelSD.h:104
float phiAtEntry_
Definition: PPSPixelSD.h:101
float tof_
Definition: PPSPixelSD.h:98
float vx_
Definition: PPSPixelSD.h:102
int primaryID_
Definition: PPSPixelSD.h:109
G4Track * theTrack_
Definition: PPSPixelSD.h:83
float vy_
Definition: PPSPixelSD.h:103
float thetaAtEntry_
Definition: PPSPixelSD.h:100
G4ThreeVector theLocalExitPoint_
Definition: PPSPixelSD.h:92
float incidentEnergy_
Definition: PPSPixelSD.h:96
float pabs_
Definition: PPSPixelSD.h:97
double tSlice_
Definition: PPSPixelSD.h:94
int tSliceID_
Definition: PPSPixelSD.h:111
G4ThreeVector hitPoint_
Definition: PPSPixelSD.h:89
G4StepPoint * postStepPoint_
Definition: PPSPixelSD.h:88
short particleType_
Definition: PPSPixelSD.h:113
uint32_t unitID_
Definition: PPSPixelSD.h:106
G4ThreeVector theLocalEntryPoint_
Definition: PPSPixelSD.h:91
#define LogDebug(id)

◆ storeHit()

void PPSPixelSD::storeHit ( PPSPixelG4Hit hit)
private

Definition at line 260 of file PPSPixelSD.cc.

References theHC_.

Referenced by createNewHit().

260  {
261  if (hit == nullptr) {
262  edm::LogWarning("PPSSim") << "PPSPixelSD: hit to be stored is NULL !!";
263  return;
264  }
265 
266  theHC_->insert(hit);
267 }
PPSPixelG4HitCollection * theHC_
Definition: PPSPixelSD.h:80
Log< level::Warning, false > LogWarning

◆ update() [1/2]

void PPSPixelSD::update ( const BeginOfEvent )
overrideprotectedvirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfEvent *>.

Definition at line 129 of file PPSPixelSD.cc.

References clearHits(), eventno_, and LogDebug.

Referenced by progressbar.ProgressBar::__next__(), MatrixUtil.Matrix::__setitem__(), MatrixUtil.Steps::__setitem__(), progressbar.ProgressBar::finish(), and MatrixUtil.Steps::overwrite().

129  {
130  LogDebug("PPSSim") << "PPSPixelSD: dispatched BeginOfEvent for " << GetName();
131  clearHits();
132  eventno_ = (*i)()->GetEventID();
133 }
void clearHits() override
Definition: PPSPixelSD.cc:137
int eventno_
Definition: PPSPixelSD.h:112
#define LogDebug(id)

◆ update() [2/2]

void PPSPixelSD::update ( const ::EndOfEvent )
overrideprotected

◆ updateHit()

void PPSPixelSD::updateHit ( )
private

Definition at line 249 of file PPSPixelSD.cc.

References currentHit_, eloss_, LogDebug, postStepPoint_, previousUnitID_, PPSPixelG4Hit::setEnergyLoss(), tsID_, tSliceID_, and unitID_.

Referenced by createNewHit(), and hitExists().

249  {
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 }
uint32_t previousUnitID_
Definition: PPSPixelSD.h:107
PPSPixelG4Hit * currentHit_
Definition: PPSPixelSD.h:82
int tSliceID_
Definition: PPSPixelSD.h:111
G4StepPoint * postStepPoint_
Definition: PPSPixelSD.h:88
uint32_t unitID_
Definition: PPSPixelSD.h:106
void setEnergyLoss(float e)
double eloss_
Definition: PPSPixelSD.h:95
#define LogDebug(id)

Member Data Documentation

◆ currentHit_

PPSPixelG4Hit* PPSPixelSD::currentHit_ = nullptr
private

Definition at line 82 of file PPSPixelSD.h.

Referenced by createNewHit(), hitExists(), and updateHit().

◆ currentPV_

G4VPhysicalVolume* PPSPixelSD::currentPV_ = nullptr
private

Definition at line 84 of file PPSPixelSD.h.

Referenced by createNewHit(), and stepInfo().

◆ eloss_

double PPSPixelSD::eloss_ = 0.0
private

Definition at line 95 of file PPSPixelSD.h.

Referenced by ProcessHits(), and updateHit().

◆ eventno_

int PPSPixelSD::eventno_ = 0
private

Definition at line 112 of file PPSPixelSD.h.

Referenced by update().

◆ exitPoint_

G4ThreeVector PPSPixelSD::exitPoint_
private

Definition at line 90 of file PPSPixelSD.h.

Referenced by createNewHit(), and stepInfo().

◆ hcID_

G4int PPSPixelSD::hcID_ = -1
private

Definition at line 85 of file PPSPixelSD.h.

Referenced by Initialize().

◆ hitPoint_

G4ThreeVector PPSPixelSD::hitPoint_
private

Definition at line 89 of file PPSPixelSD.h.

Referenced by createNewHit(), and stepInfo().

◆ incidentEnergy_

float PPSPixelSD::incidentEnergy_ = 0.f
private

Definition at line 96 of file PPSPixelSD.h.

Referenced by createNewHit(), and stepInfo().

◆ numberingScheme_

std::unique_ptr<PPSVDetectorOrganization> PPSPixelSD::numberingScheme_
private

Definition at line 78 of file PPSPixelSD.h.

Referenced by PPSPixelSD(), and setDetUnitId().

◆ pabs_

float PPSPixelSD::pabs_ = 0.f
private

Definition at line 97 of file PPSPixelSD.h.

Referenced by createNewHit(), and stepInfo().

◆ parentID_

int PPSPixelSD::parentID_ = 0
private

Definition at line 110 of file PPSPixelSD.h.

Referenced by createNewHit(), and stepInfo().

◆ particleType_

short PPSPixelSD::particleType_ = 0
private

Definition at line 113 of file PPSPixelSD.h.

Referenced by createNewHit(), and stepInfo().

◆ phiAtEntry_

float PPSPixelSD::phiAtEntry_ = 0.f
private

Definition at line 101 of file PPSPixelSD.h.

Referenced by createNewHit(), and stepInfo().

◆ postStepPoint_

G4StepPoint* PPSPixelSD::postStepPoint_ = nullptr
private

Definition at line 88 of file PPSPixelSD.h.

Referenced by stepInfo(), and updateHit().

◆ preStepPoint_

G4StepPoint* PPSPixelSD::preStepPoint_ = nullptr
private

Definition at line 87 of file PPSPixelSD.h.

Referenced by setToLocal(), and stepInfo().

◆ previousUnitID_

uint32_t PPSPixelSD::previousUnitID_ = 0
private

Definition at line 107 of file PPSPixelSD.h.

Referenced by hitExists(), and updateHit().

◆ primaryID_

int PPSPixelSD::primaryID_ = 0
private

Definition at line 109 of file PPSPixelSD.h.

Referenced by createNewHit(), and stepInfo().

◆ slave_

std::unique_ptr<TrackingSlaveSD> PPSPixelSD::slave_
private

Definition at line 77 of file PPSPixelSD.h.

Referenced by clearHits(), EndOfEvent(), fillHits(), and PPSPixelSD().

◆ theHC_

PPSPixelG4HitCollection* PPSPixelSD::theHC_ = nullptr
private

Definition at line 80 of file PPSPixelSD.h.

Referenced by EndOfEvent(), hitExists(), Initialize(), PrintAll(), and storeHit().

◆ theLocalEntryPoint_

G4ThreeVector PPSPixelSD::theLocalEntryPoint_
private

Definition at line 91 of file PPSPixelSD.h.

Referenced by createNewHit(), and stepInfo().

◆ theLocalExitPoint_

G4ThreeVector PPSPixelSD::theLocalExitPoint_
private

Definition at line 92 of file PPSPixelSD.h.

Referenced by createNewHit(), and stepInfo().

◆ theManager_

const SimTrackManager* PPSPixelSD::theManager_
private

Definition at line 81 of file PPSPixelSD.h.

◆ thetaAtEntry_

float PPSPixelSD::thetaAtEntry_ = 0.f
private

Definition at line 100 of file PPSPixelSD.h.

Referenced by createNewHit(), and stepInfo().

◆ theTrack_

G4Track* PPSPixelSD::theTrack_ = nullptr
private

Definition at line 83 of file PPSPixelSD.h.

Referenced by createNewHit(), ProcessHits(), and stepInfo().

◆ tof_

float PPSPixelSD::tof_ = 0.f
private

Definition at line 98 of file PPSPixelSD.h.

Referenced by createNewHit(), and stepInfo().

◆ tsID_

int PPSPixelSD::tsID_ = 0
private

Definition at line 108 of file PPSPixelSD.h.

Referenced by hitExists(), Initialize(), and updateHit().

◆ tSlice_

double PPSPixelSD::tSlice_ = 0.0
private

Definition at line 94 of file PPSPixelSD.h.

Referenced by createNewHit(), and stepInfo().

◆ tSliceID_

int PPSPixelSD::tSliceID_ = 0
private

Definition at line 111 of file PPSPixelSD.h.

Referenced by createNewHit(), hitExists(), stepInfo(), and updateHit().

◆ unitID_

uint32_t PPSPixelSD::unitID_ = 0
private

Definition at line 106 of file PPSPixelSD.h.

Referenced by createNewHit(), hitExists(), stepInfo(), and updateHit().

◆ vx_

float PPSPixelSD::vx_ = 0.f
private

Definition at line 102 of file PPSPixelSD.h.

Referenced by createNewHit(), and stepInfo().

◆ vy_

float PPSPixelSD::vy_ = 0.f
private

Definition at line 103 of file PPSPixelSD.h.

Referenced by createNewHit(), and stepInfo().

◆ vz_

float PPSPixelSD::vz_ = 0.f
private

Definition at line 104 of file PPSPixelSD.h.

Referenced by createNewHit(), and stepInfo().