CMS 3D CMS Logo

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

#include <PPSDiamondSD.h>

Inheritance diagram for PPSDiamondSD:
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
 
 PPSDiamondSD (const std::string &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
 
void PrintAll () override
 
bool ProcessHits (G4Step *step, G4TouchableHistory *tHistory) override
 
uint32_t setDetUnitId (const G4Step *) override
 
 ~PPSDiamondSD () 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 importInfoToHit ()
 
void printHitInfo ()
 
G4ThreeVector setToLocal (const G4ThreeVector &globalPoint)
 
void stepInfo (const G4Step *aStep)
 
void storeHit (PPSDiamondG4Hit *)
 

Private Attributes

PPSDiamondG4HitcurrentHit_ = nullptr
 
G4VPhysicalVolume * currentPV_ = nullptr
 
double eloss_ = 0.0
 
int eventno_ = 0
 
G4ThreeVector exitPoint_
 
double Globaltimehit_
 
G4int hcID_ = -1
 
G4ThreeVector hitPoint_
 
double incidentEnergy_ = 0.0
 
std::unique_ptr< PPSVDetectorOrganizationnumberingScheme_
 
double pabs_ = 0.0
 
G4int parentID_ = 0
 
short particleType_ = 0
 
double phiAtEntry_ = 0.0
 
G4StepPoint * postStepPoint_ = nullptr
 
G4StepPoint * preStepPoint_ = nullptr
 
G4int primaryID_ = 0
 
std::unique_ptr< TrackingSlaveSDslave_
 
double theglobaltimehit_
 
PPSDiamondG4HitCollectiontheHC_ = nullptr
 
G4ThreeVector theLocalEntryPoint_
 
G4ThreeVector theLocalExitPoint_
 
double thePx_ = 0.0
 
double thePy_ = 0.0
 
double thePz_ = 0.0
 
double thetaAtEntry_ = 0.0
 
G4Track * theTrack_ = nullptr
 
double tof_ = 0.0
 
G4double tSlice_ = 0.0
 
G4int tSliceID_ = 0
 
unsigned int unitID_ = 0
 
int verbosity_
 
double vx_ = 0.0
 
double vy_ = 0.0
 
double vz_ = 0.0
 

Additional Inherited Members

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

Detailed Description

Definition at line 29 of file PPSDiamondSD.h.

Constructor & Destructor Documentation

◆ PPSDiamondSD()

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

Definition at line 29 of file PPSDiamondSD.cc.

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

33  : SensitiveTkDetector(pname, clg) {
34  collectionName.insert(pname);
35  edm::ParameterSet m_Anal = p.getParameter<edm::ParameterSet>("PPSDiamondSD");
36  verbosity_ = m_Anal.getParameter<int>("Verbosity");
37 
38  LogDebug("PPSSimDiamond") << "*******************************************************\n"
39  << "* *\n"
40  << "* Constructing a PPSDiamondSD with name " << pname << "\n"
41  << "* *\n"
42  << "*******************************************************"
43  << "\n";
44 
45  slave_ = std::make_unique<TrackingSlaveSD>(pname);
46 
47  if (pname == "CTPPSTimingHits") {
48  numberingScheme_ = std::make_unique<PPSDiamondOrganization>();
49  edm::LogVerbatim("PPSSimDiamond") << "Find CTPPSTimingHits as name";
50  } else {
51  edm::LogError("PPSSimDiamond") << "PPSDiamondSD: ReadoutName " << pname << " not supported";
52  }
53 
54  edm::LogVerbatim("PPSSimDiamond") << "PPSDiamondSD: Instantiation completed for " << pname;
55 }
Log< level::Info, true > LogVerbatim
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
std::unique_ptr< TrackingSlaveSD > slave_
Definition: PPSDiamondSD.h:59
std::unique_ptr< PPSVDetectorOrganization > numberingScheme_
Definition: PPSDiamondSD.h:60
Log< level::Error, false > LogError
SensitiveTkDetector(const std::string &iname, const SensitiveDetectorCatalog &clg)
#define LogDebug(id)

◆ ~PPSDiamondSD()

PPSDiamondSD::~PPSDiamondSD ( )
override

Definition at line 57 of file PPSDiamondSD.cc.

57 {}

Member Function Documentation

◆ clearHits()

void PPSDiamondSD::clearHits ( )
overridevirtual

Implements SensitiveDetector.

Definition at line 244 of file PPSDiamondSD.cc.

References slave_.

Referenced by update().

244 { slave_->Initialize(); }
std::unique_ptr< TrackingSlaveSD > slave_
Definition: PPSDiamondSD.h:59

◆ EndOfEvent()

void PPSDiamondSD::EndOfEvent ( G4HCofThisEvent *  eventHC)
override

Definition at line 205 of file PPSDiamondSD.cc.

References PPSDiamondG4Hit::energyLoss(), mps_splice::entry, beamvalidation::exit(), createfilelist::int, dqmiolumiharvest::j, PPSDiamondG4Hit::localEntry(), PPSDiamondG4Hit::localExit(), PPSDiamondG4Hit::p(), PPSDiamondG4Hit::particleType(), PPSDiamondG4Hit::phiAtEntry(), slave_, theHC_, PPSDiamondG4Hit::thetaAtEntry(), PPSDiamondG4Hit::tof(), PPSDiamondG4Hit::trackID(), and PPSDiamondG4Hit::unitID().

205  {
206  // here we loop over transient hits and make them persistent
207  for (unsigned int j = 0; j < (unsigned int)theHC_->entries(); ++j) {
208  PPSDiamondG4Hit* aHit = (*theHC_)[j];
209 
210  Local3DPoint entry(aHit->localEntry().x(), aHit->localEntry().y(), aHit->localEntry().z());
211  Local3DPoint exit(aHit->localExit().x(), aHit->localExit().y(), aHit->localExit().z());
212  slave_->processHits(PSimHit(entry,
213  exit,
214  aHit->p(),
215  aHit->tof(),
216  aHit->energyLoss(),
217  aHit->particleType(),
218  aHit->unitID(),
219  aHit->trackID(),
220  aHit->thetaAtEntry(),
221  aHit->phiAtEntry()));
222  }
223 }
double thetaAtEntry() const
std::unique_ptr< TrackingSlaveSD > slave_
Definition: PPSDiamondSD.h:59
double phiAtEntry() const
double energyLoss() const
double p() const
int unitID() const
double tof() const
const G4ThreeVector & localExit() const
PPSDiamondG4HitCollection * theHC_
Definition: PPSDiamondSD.h:62
const G4ThreeVector & localEntry() const
unsigned int trackID() const
int particleType() const
def exit(msg="")

◆ fillHits()

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

Implements SensitiveTkDetector.

Definition at line 230 of file PPSDiamondSD.cc.

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

230  {
231  if (slave_->name() == n)
232  c = slave_->hits();
233 }
std::unique_ptr< TrackingSlaveSD > slave_
Definition: PPSDiamondSD.h:59

◆ importInfoToHit()

void PPSDiamondSD::importInfoToHit ( )
private

Definition at line 168 of file PPSDiamondSD.cc.

References currentHit_, eloss_, exitPoint_, hitPoint_, incidentEnergy_, LogDebug, pabs_, parentID_, particleType_, phiAtEntry_, primaryID_, PPSDiamondG4Hit::setEnergyLoss(), PPSDiamondG4Hit::setEntry(), PPSDiamondG4Hit::setExit(), PPSDiamondG4Hit::setGlobalTimehit(), PPSDiamondG4Hit::setIncidentEnergy(), PPSDiamondG4Hit::setLocalEntry(), PPSDiamondG4Hit::setLocalExit(), PPSDiamondG4Hit::setP(), PPSDiamondG4Hit::setParentId(), PPSDiamondG4Hit::setParticleType(), PPSDiamondG4Hit::setPhiAtEntry(), PPSDiamondG4Hit::setPx(), PPSDiamondG4Hit::setPy(), PPSDiamondG4Hit::setPz(), PPSDiamondG4Hit::setThetaAtEntry(), PPSDiamondG4Hit::setTimeSlice(), PPSDiamondG4Hit::setTof(), PPSDiamondG4Hit::setTrackID(), PPSDiamondG4Hit::setUnitID(), PPSDiamondG4Hit::setVx(), PPSDiamondG4Hit::setVy(), PPSDiamondG4Hit::setVz(), storeHit(), theLocalEntryPoint_, theLocalExitPoint_, thePx_, thePy_, thePz_, thetaAtEntry_, tof_, tSlice_, unitID_, vx_, vy_, and vz_.

Referenced by ProcessHits().

168  {
192 
194  LogDebug("PPSSimDiamond") << "STORED HIT IN: " << unitID_ << "\n";
195 }
unsigned int unitID_
Definition: PPSDiamondSD.h:97
void setP(double e)
PPSDiamondG4Hit * currentHit_
Definition: PPSDiamondSD.h:63
void setTof(double e)
G4ThreeVector hitPoint_
Definition: PPSDiamondSD.h:74
void setVy(double p)
void setPz(double p)
G4int primaryID_
Definition: PPSDiamondSD.h:67
void setEntry(const G4ThreeVector &xyz)
G4ThreeVector theLocalEntryPoint_
Definition: PPSDiamondSD.h:76
void setIncidentEnergy(double e)
void setThetaAtEntry(double t)
void setVz(double p)
G4ThreeVector theLocalExitPoint_
Definition: PPSDiamondSD.h:77
void setPy(double p)
double eloss_
Definition: PPSDiamondSD.h:85
void setUnitID(unsigned int i)
void setEnergyLoss(double e)
void setExit(const G4ThreeVector &xyz)
void setVx(double p)
G4ThreeVector exitPoint_
Definition: PPSDiamondSD.h:75
double thePy_
Definition: PPSDiamondSD.h:82
G4double tSlice_
Definition: PPSDiamondSD.h:70
void storeHit(PPSDiamondG4Hit *)
short particleType_
Definition: PPSDiamondSD.h:100
void setPhiAtEntry(double f)
void setTrackID(int i)
double pabs_
Definition: PPSDiamondSD.h:80
void setParticleType(short i)
void setParentId(int p)
void setGlobalTimehit(double h)
void setLocalEntry(const G4ThreeVector &theLocalEntryPoint)
double thetaAtEntry_
Definition: PPSDiamondSD.h:87
double thePz_
Definition: PPSDiamondSD.h:83
double phiAtEntry_
Definition: PPSDiamondSD.h:88
G4int parentID_
Definition: PPSDiamondSD.h:68
double incidentEnergy_
Definition: PPSDiamondSD.h:79
double thePx_
Definition: PPSDiamondSD.h:81
void setLocalExit(const G4ThreeVector &theLocalExitPoint)
void setPx(double p)
void setTimeSlice(double d)
#define LogDebug(id)

◆ Initialize()

void PPSDiamondSD::Initialize ( G4HCofThisEvent *  HCE)
override

Definition at line 59 of file PPSDiamondSD.cc.

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

59  {
60  LogDebug("PPSSimDiamond") << "PPSDiamondSD : Initialize called for " << GetName();
61 
62  theHC_ = new PPSDiamondG4HitCollection(GetName(), collectionName[0]);
63  G4SDManager::GetSDMpointer()->AddNewCollection(GetName(), collectionName[0]);
64  if (hcID_ < 0)
65  hcID_ = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]);
66  HCE->AddHitsCollection(hcID_, theHC_);
67 }
PPSDiamondG4HitCollection * theHC_
Definition: PPSDiamondSD.h:62
G4THitsCollection< PPSDiamondG4Hit > PPSDiamondG4HitCollection
#define LogDebug(id)

◆ PrintAll()

void PPSDiamondSD::PrintAll ( )
override

Definition at line 225 of file PPSDiamondSD.cc.

References LogDebug, and theHC_.

225  {
226  LogDebug("PPSSimDiamond") << "PPSDiamond: Collection " << theHC_->GetName() << "\n";
227  theHC_->PrintAllHits();
228 }
PPSDiamondG4HitCollection * theHC_
Definition: PPSDiamondSD.h:62
#define LogDebug(id)

◆ printHitInfo()

void PPSDiamondSD::printHitInfo ( )
private

Definition at line 69 of file PPSDiamondSD.cc.

References currentPV_, eloss_, exitPoint_, hitPoint_, LogDebug, parentID_, postStepPoint_, primaryID_, theTrack_, tSliceID_, and unitID_.

69  {
70  LogDebug("PPSSimDiamond") << theTrack_->GetDefinition()->GetParticleName() << " PPS_Timing_SD CreateNewHit for"
71  << " PV " << currentPV_->GetName() << " PVid = " << currentPV_->GetCopyNo() << " Unit "
72  << unitID_ << "\n";
73  LogDebug("PPSSimDiamond") << " primary " << primaryID_ << " time slice " << tSliceID_ << " of energy "
74  << theTrack_->GetTotalEnergy() << " Eloss=" << eloss_ << " position pre: " << hitPoint_
75  << " post: " << exitPoint_;
76  LogDebug("PPSSimDiamond") << " positions "
77  << "(" << postStepPoint_->GetPosition().x() << "," << postStepPoint_->GetPosition().y()
78  << "," << postStepPoint_->GetPosition().z() << ")"
79  << " For Track " << theTrack_->GetTrackID() << " which is a "
80  << theTrack_->GetDefinition()->GetParticleName() << " ParentID is "
81  << theTrack_->GetParentID() << "\n";
82 
83  if (parentID_ == 0) {
84  LogDebug("PPSSimDiamond") << " primary particle ";
85  } else {
86  LogDebug("PPSSimDiamond") << " daughter of part. " << parentID_;
87  }
88 
89  LogDebug("PPSSimDiamond") << " and created by ";
90 
91  if (theTrack_->GetCreatorProcess() != nullptr)
92  LogDebug("PPSSimDiamond") << theTrack_->GetCreatorProcess()->GetProcessName();
93  else
94  LogDebug("PPSSimDiamond") << "NO process";
95 }
G4VPhysicalVolume * currentPV_
Definition: PPSDiamondSD.h:65
unsigned int unitID_
Definition: PPSDiamondSD.h:97
G4ThreeVector hitPoint_
Definition: PPSDiamondSD.h:74
G4int primaryID_
Definition: PPSDiamondSD.h:67
double eloss_
Definition: PPSDiamondSD.h:85
G4ThreeVector exitPoint_
Definition: PPSDiamondSD.h:75
G4Track * theTrack_
Definition: PPSDiamondSD.h:64
G4StepPoint * postStepPoint_
Definition: PPSDiamondSD.h:73
G4int parentID_
Definition: PPSDiamondSD.h:68
G4int tSliceID_
Definition: PPSDiamondSD.h:69
#define LogDebug(id)

◆ ProcessHits()

bool PPSDiamondSD::ProcessHits ( G4Step *  step,
G4TouchableHistory *  tHistory 
)
overridevirtual

Implements SensitiveDetector.

Definition at line 97 of file PPSDiamondSD.cc.

References eloss_, importInfoToHit(), LogDebug, stepInfo(), and theTrack_.

97  {
98  eloss_ = aStep->GetTotalEnergyDeposit();
99  if (eloss_ > 0.0) {
100  eloss_ /= GeV;
101  stepInfo(aStep);
102  if (theTrack_->GetDefinition()->GetPDGEncoding() == 2212) {
103  LogDebug("PPSSimDiamond") << "PPSSimDiamond : eloss=" << eloss_;
104  importInfoToHit(); //in addition to import info to hit it STORE hit as well
105  LogDebug("PPSSimDiamond") << " information imported to the hit ";
106  }
107  }
108  return true;
109 }
void stepInfo(const G4Step *aStep)
double eloss_
Definition: PPSDiamondSD.h:85
G4Track * theTrack_
Definition: PPSDiamondSD.h:64
void importInfoToHit()
#define LogDebug(id)

◆ setDetUnitId()

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

Implements SensitiveDetector.

Definition at line 154 of file PPSDiamondSD.cc.

References numberingScheme_.

Referenced by stepInfo().

154  {
155  return (numberingScheme_ == nullptr ? 0 : numberingScheme_->unitID(aStep));
156 }
std::unique_ptr< PPSVDetectorOrganization > numberingScheme_
Definition: PPSDiamondSD.h:60

◆ setToLocal()

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

Definition at line 197 of file PPSDiamondSD.cc.

References preStepPoint_.

Referenced by stepInfo().

197  {
198  G4ThreeVector localPoint;
199  const G4VTouchable* touch = preStepPoint_->GetTouchable();
200  localPoint = touch->GetHistory()->GetTopTransform().TransformPoint(global);
201 
202  return localPoint;
203 }
G4StepPoint * preStepPoint_
Definition: PPSDiamondSD.h:72

◆ stepInfo()

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

Definition at line 111 of file PPSDiamondSD.cc.

References SensitiveDetector::ConvertToLocal3DPoint(), currentPV_, exitPoint_, hitPoint_, incidentEnergy_, createfilelist::int, LogDebug, pabs_, parentID_, particleType_, PV3DBase< T, PVType, FrameType >::phi(), phiAtEntry_, postStepPoint_, preStepPoint_, primaryID_, setDetUnitId(), setToLocal(), theLocalEntryPoint_, theLocalExitPoint_, thePx_, thePy_, thePz_, PV3DBase< T, PVType, FrameType >::theta(), thetaAtEntry_, theTrack_, tof_, tSlice_, tSliceID_, unitID_, verbosity_, vx_, vy_, and vz_.

Referenced by ProcessHits().

111  {
112  theTrack_ = aStep->GetTrack();
113  preStepPoint_ = aStep->GetPreStepPoint();
114  postStepPoint_ = aStep->GetPostStepPoint();
115  hitPoint_ = preStepPoint_->GetPosition();
116  exitPoint_ = postStepPoint_->GetPosition();
117  currentPV_ = preStepPoint_->GetPhysicalVolume();
120  tof_ = preStepPoint_->GetGlobalTime() / nanosecond;
121  incidentEnergy_ = preStepPoint_->GetTotalEnergy() / eV;
122  tSlice_ = postStepPoint_->GetGlobalTime() / nanosecond;
123  tSliceID_ = (int)tSlice_;
124  unitID_ = setDetUnitId(aStep);
125 
126  if (verbosity_)
127  LogDebug("PPSSimDiamond") << "UNIT " << unitID_;
128 
129  primaryID_ = theTrack_->GetTrackID();
130  parentID_ = theTrack_->GetParentID();
131 
132  incidentEnergy_ = theTrack_->GetTotalEnergy() / GeV;
133 
134  pabs_ = preStepPoint_->GetMomentum().mag() / GeV;
135  thePx_ = preStepPoint_->GetMomentum().x() / GeV;
136  thePy_ = preStepPoint_->GetMomentum().y() / GeV;
137  thePz_ = preStepPoint_->GetMomentum().z() / GeV;
138  particleType_ = theTrack_->GetDefinition()->GetPDGEncoding();
139 
140  //corrected phi and theta treatment
141  G4ThreeVector gmd = preStepPoint_->GetMomentumDirection();
142  // convert it to local frame
143  G4ThreeVector lmd =
144  ((G4TouchableHistory*)(preStepPoint_->GetTouchable()))->GetHistory()->GetTopTransform().TransformAxis(gmd);
146  thetaAtEntry_ = lnmd.theta();
147  phiAtEntry_ = lnmd.phi();
148 
149  vx_ = theTrack_->GetVertexPosition().x() / mm;
150  vy_ = theTrack_->GetVertexPosition().y() / mm;
151  vz_ = theTrack_->GetVertexPosition().z() / mm;
152 }
G4VPhysicalVolume * currentPV_
Definition: PPSDiamondSD.h:65
unsigned int unitID_
Definition: PPSDiamondSD.h:97
G4ThreeVector setToLocal(const G4ThreeVector &globalPoint)
G4ThreeVector hitPoint_
Definition: PPSDiamondSD.h:74
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
G4int primaryID_
Definition: PPSDiamondSD.h:67
G4ThreeVector theLocalEntryPoint_
Definition: PPSDiamondSD.h:76
G4ThreeVector theLocalExitPoint_
Definition: PPSDiamondSD.h:77
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
G4ThreeVector exitPoint_
Definition: PPSDiamondSD.h:75
double thePy_
Definition: PPSDiamondSD.h:82
G4double tSlice_
Definition: PPSDiamondSD.h:70
short particleType_
Definition: PPSDiamondSD.h:100
uint32_t setDetUnitId(const G4Step *) override
G4Track * theTrack_
Definition: PPSDiamondSD.h:64
double pabs_
Definition: PPSDiamondSD.h:80
double thetaAtEntry_
Definition: PPSDiamondSD.h:87
double thePz_
Definition: PPSDiamondSD.h:83
G4StepPoint * postStepPoint_
Definition: PPSDiamondSD.h:73
double phiAtEntry_
Definition: PPSDiamondSD.h:88
G4int parentID_
Definition: PPSDiamondSD.h:68
double incidentEnergy_
Definition: PPSDiamondSD.h:79
G4StepPoint * preStepPoint_
Definition: PPSDiamondSD.h:72
double thePx_
Definition: PPSDiamondSD.h:81
G4int tSliceID_
Definition: PPSDiamondSD.h:69
Geom::Theta< T > theta() const
Definition: PV3DBase.h:72
#define LogDebug(id)

◆ storeHit()

void PPSDiamondSD::storeHit ( PPSDiamondG4Hit hit)
private

Definition at line 158 of file PPSDiamondSD.cc.

References LogDebug, theHC_, and verbosity_.

Referenced by importInfoToHit().

158  {
159  if (hit == nullptr) {
160  if (verbosity_)
161  LogDebug("PPSSimDiamond") << "PPSDiamond: hit to be stored is NULL !!";
162  return;
163  }
164 
165  theHC_->insert(hit);
166 }
PPSDiamondG4HitCollection * theHC_
Definition: PPSDiamondSD.h:62
#define LogDebug(id)

◆ update() [1/2]

void PPSDiamondSD::update ( const BeginOfEvent )
overrideprotectedvirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfEvent *>.

Definition at line 235 of file PPSDiamondSD.cc.

References clearHits(), eventno_, and LogDebug.

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

235  {
236  LogDebug("PPSSimDiamond") << " Dispatched BeginOfEvent !"
237  << "\n";
238  clearHits();
239  eventno_ = (*i)()->GetEventID();
240 }
void clearHits() override
#define LogDebug(id)

◆ update() [2/2]

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

Member Data Documentation

◆ currentHit_

PPSDiamondG4Hit* PPSDiamondSD::currentHit_ = nullptr
private

Definition at line 63 of file PPSDiamondSD.h.

Referenced by importInfoToHit().

◆ currentPV_

G4VPhysicalVolume* PPSDiamondSD::currentPV_ = nullptr
private

Definition at line 65 of file PPSDiamondSD.h.

Referenced by printHitInfo(), and stepInfo().

◆ eloss_

double PPSDiamondSD::eloss_ = 0.0
private

Definition at line 85 of file PPSDiamondSD.h.

Referenced by importInfoToHit(), printHitInfo(), and ProcessHits().

◆ eventno_

int PPSDiamondSD::eventno_ = 0
private

Definition at line 99 of file PPSDiamondSD.h.

Referenced by update().

◆ exitPoint_

G4ThreeVector PPSDiamondSD::exitPoint_
private

Definition at line 75 of file PPSDiamondSD.h.

Referenced by importInfoToHit(), printHitInfo(), and stepInfo().

◆ Globaltimehit_

double PPSDiamondSD::Globaltimehit_
private

Definition at line 94 of file PPSDiamondSD.h.

◆ hcID_

G4int PPSDiamondSD::hcID_ = -1
private

Definition at line 66 of file PPSDiamondSD.h.

Referenced by Initialize().

◆ hitPoint_

G4ThreeVector PPSDiamondSD::hitPoint_
private

Definition at line 74 of file PPSDiamondSD.h.

Referenced by importInfoToHit(), printHitInfo(), and stepInfo().

◆ incidentEnergy_

double PPSDiamondSD::incidentEnergy_ = 0.0
private

Definition at line 79 of file PPSDiamondSD.h.

Referenced by importInfoToHit(), and stepInfo().

◆ numberingScheme_

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

Definition at line 60 of file PPSDiamondSD.h.

Referenced by PPSDiamondSD(), and setDetUnitId().

◆ pabs_

double PPSDiamondSD::pabs_ = 0.0
private

Definition at line 80 of file PPSDiamondSD.h.

Referenced by importInfoToHit(), and stepInfo().

◆ parentID_

G4int PPSDiamondSD::parentID_ = 0
private

Definition at line 68 of file PPSDiamondSD.h.

Referenced by importInfoToHit(), printHitInfo(), and stepInfo().

◆ particleType_

short PPSDiamondSD::particleType_ = 0
private

Definition at line 100 of file PPSDiamondSD.h.

Referenced by importInfoToHit(), and stepInfo().

◆ phiAtEntry_

double PPSDiamondSD::phiAtEntry_ = 0.0
private

Definition at line 88 of file PPSDiamondSD.h.

Referenced by importInfoToHit(), and stepInfo().

◆ postStepPoint_

G4StepPoint* PPSDiamondSD::postStepPoint_ = nullptr
private

Definition at line 73 of file PPSDiamondSD.h.

Referenced by printHitInfo(), and stepInfo().

◆ preStepPoint_

G4StepPoint* PPSDiamondSD::preStepPoint_ = nullptr
private

Definition at line 72 of file PPSDiamondSD.h.

Referenced by setToLocal(), and stepInfo().

◆ primaryID_

G4int PPSDiamondSD::primaryID_ = 0
private

Definition at line 67 of file PPSDiamondSD.h.

Referenced by importInfoToHit(), printHitInfo(), and stepInfo().

◆ slave_

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

Definition at line 59 of file PPSDiamondSD.h.

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

◆ theglobaltimehit_

double PPSDiamondSD::theglobaltimehit_
private

Definition at line 95 of file PPSDiamondSD.h.

◆ theHC_

PPSDiamondG4HitCollection* PPSDiamondSD::theHC_ = nullptr
private

Definition at line 62 of file PPSDiamondSD.h.

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

◆ theLocalEntryPoint_

G4ThreeVector PPSDiamondSD::theLocalEntryPoint_
private

Definition at line 76 of file PPSDiamondSD.h.

Referenced by importInfoToHit(), and stepInfo().

◆ theLocalExitPoint_

G4ThreeVector PPSDiamondSD::theLocalExitPoint_
private

Definition at line 77 of file PPSDiamondSD.h.

Referenced by importInfoToHit(), and stepInfo().

◆ thePx_

double PPSDiamondSD::thePx_ = 0.0
private

Definition at line 81 of file PPSDiamondSD.h.

Referenced by importInfoToHit(), and stepInfo().

◆ thePy_

double PPSDiamondSD::thePy_ = 0.0
private

Definition at line 82 of file PPSDiamondSD.h.

Referenced by importInfoToHit(), and stepInfo().

◆ thePz_

double PPSDiamondSD::thePz_ = 0.0
private

Definition at line 83 of file PPSDiamondSD.h.

Referenced by importInfoToHit(), and stepInfo().

◆ thetaAtEntry_

double PPSDiamondSD::thetaAtEntry_ = 0.0
private

Definition at line 87 of file PPSDiamondSD.h.

Referenced by importInfoToHit(), and stepInfo().

◆ theTrack_

G4Track* PPSDiamondSD::theTrack_ = nullptr
private

Definition at line 64 of file PPSDiamondSD.h.

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

◆ tof_

double PPSDiamondSD::tof_ = 0.0
private

Definition at line 84 of file PPSDiamondSD.h.

Referenced by importInfoToHit(), and stepInfo().

◆ tSlice_

G4double PPSDiamondSD::tSlice_ = 0.0
private

Definition at line 70 of file PPSDiamondSD.h.

Referenced by importInfoToHit(), and stepInfo().

◆ tSliceID_

G4int PPSDiamondSD::tSliceID_ = 0
private

Definition at line 69 of file PPSDiamondSD.h.

Referenced by printHitInfo(), and stepInfo().

◆ unitID_

unsigned int PPSDiamondSD::unitID_ = 0
private

Definition at line 97 of file PPSDiamondSD.h.

Referenced by importInfoToHit(), printHitInfo(), and stepInfo().

◆ verbosity_

int PPSDiamondSD::verbosity_
private

Definition at line 98 of file PPSDiamondSD.h.

Referenced by PPSDiamondSD(), stepInfo(), and storeHit().

◆ vx_

double PPSDiamondSD::vx_ = 0.0
private

Definition at line 90 of file PPSDiamondSD.h.

Referenced by importInfoToHit(), and stepInfo().

◆ vy_

double PPSDiamondSD::vy_ = 0.0
private

Definition at line 91 of file PPSDiamondSD.h.

Referenced by importInfoToHit(), and stepInfo().

◆ vz_

double PPSDiamondSD::vz_ = 0.0
private

Definition at line 92 of file PPSDiamondSD.h.

Referenced by importInfoToHit(), and stepInfo().