CMS 3D CMS Logo

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

#include <TkAccumulatingSensitiveDetector.h>

Inheritance diagram for TkAccumulatingSensitiveDetector:
SensitiveTkDetector Observer< const BeginOfEvent *> Observer< const BeginOfTrack *> Observer< const BeginOfJob *> SensitiveDetector

Public Member Functions

void clearHits () override
 
void EndOfEvent (G4HCofThisEvent *) override
 
void fillHits (edm::PSimHitContainer &, const std::string &) override
 
bool ProcessHits (G4Step *, G4TouchableHistory *) override
 
uint32_t setDetUnitId (const G4Step *) override
 
 TkAccumulatingSensitiveDetector (const std::string &, const GeometricDet *, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
 
 ~TkAccumulatingSensitiveDetector () 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)
 
 ~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 BeginOfTrack *>
 Observer ()
 
void slotForUpdate (const BeginOfTrack * iT)
 
virtual ~Observer ()
 
- Public Member Functions inherited from Observer< const BeginOfJob *>
 Observer ()
 
void slotForUpdate (const BeginOfJob * 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 BeginOfTrack *) override
 This routine will be called when the appropriate signal arrives. More...
 
void update (const BeginOfJob *) override
 This routine will be called when the appropriate signal arrives. More...
 
- 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 > &)
 

Private Member Functions

bool closeHit (const G4Step *)
 
void createHit (const G4Step *)
 
bool newHit (const G4Step *)
 
void sendHit ()
 
void updateHit (const G4Step *)
 

Private Attributes

bool allowZeroEnergyLoss
 
float energyCut
 
float energyHistoryCut
 
int eventno
 
Local3DPoint globalEntryPoint
 
Local3DPoint globalExitPoint
 
uint32_t lastId
 
int lastTrack
 
UpdatablePSimHitmySimHit
 
bool neverAccumulate
 
const G4VPhysicalVolume * oldVolume
 
const GeometricDetpDD_
 
std::string pname
 
bool printHits
 
float px
 
float py
 
float pz
 
double rTracker
 
double rTracker2
 
std::unique_ptr< TrackingSlaveSDslaveHighTof
 
std::unique_ptr< TrackingSlaveSDslaveLowTof
 
std::unique_ptr< const G4ProcessTypeEnumeratortheG4ProcTypeEnumerator
 
const SimTrackManagertheManager
 
TrackerG4SimHitNumberingSchemetheNumberingScheme
 
std::unique_ptr< FrameRotationtheRotation
 
float theTofLimit
 
double zTracker
 

Additional Inherited Members

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

Detailed Description

Definition at line 31 of file TkAccumulatingSensitiveDetector.h.

Constructor & Destructor Documentation

◆ TkAccumulatingSensitiveDetector()

TkAccumulatingSensitiveDetector::TkAccumulatingSensitiveDetector ( const std::string &  name,
const GeometricDet pDD,
const SensitiveDetectorCatalog clg,
edm::ParameterSet const &  p,
const SimTrackManager manager 
)

Definition at line 41 of file TkAccumulatingSensitiveDetector.cc.

References allowZeroEnergyLoss, energyCut, energyHistoryCut, edm::ParameterSet::getParameter(), Skims_PA_cff::name, neverAccumulate, AlCaHLTBitMon_ParallelJobs::p, printHits, rTracker, rTracker2, SensitiveDetector::setNames(), slaveHighTof, slaveLowTof, AlCaHLTBitMon_QueryRunRegistry::string, groupFilesInBlocks::temp, theG4ProcTypeEnumerator, theNumberingScheme, theRotation, theTofLimit, and zTracker.

46  : SensitiveTkDetector(name, clg),
47  pDD_(pDD),
48  theManager(manager),
49  rTracker(1200. * CLHEP::mm),
50  zTracker(3000. * CLHEP::mm),
51  mySimHit(nullptr),
52  lastId(0),
53  lastTrack(0),
54  oldVolume(nullptr),
55  px(0.0f),
56  py(0.0f),
57  pz(0.0f),
58  eventno(0),
59  pname("") {
60  edm::ParameterSet m_TrackerSD = p.getParameter<edm::ParameterSet>("TrackerSD");
61  allowZeroEnergyLoss = m_TrackerSD.getParameter<bool>("ZeroEnergyLoss");
62  neverAccumulate = m_TrackerSD.getParameter<bool>("NeverAccumulate");
63  printHits = m_TrackerSD.getParameter<bool>("PrintHits");
64  theTofLimit = m_TrackerSD.getParameter<double>("ElectronicSigmaInNanoSeconds") * 3 * CLHEP::ns; // 3 sigma
65  energyCut =
66  m_TrackerSD.getParameter<double>("EnergyThresholdForPersistencyInGeV") * CLHEP::GeV; //default must be 0.5
68  m_TrackerSD.getParameter<double>("EnergyThresholdForHistoryInGeV") * CLHEP::GeV; //default must be 0.05
70 
71  // No Rotation given in input, automagically choose one based upon the name
72  std::string rotType;
73  theRotation = std::make_unique<TrackerFrameRotation>();
74  rotType = "TrackerFrameRotation";
75 
76 #ifdef FAKEFRAMEROTATION
77  theRotation.reset(new FakeFrameRotation());
78  rotType = "FakeFrameRotation";
79 #endif
80 
81  edm::LogVerbatim("TrackerSim") << " TkAccumulatingSensitiveDetector: "
82  << " Criteria for Saving Tracker SimTracks: \n"
83  << " History: " << energyHistoryCut << " MeV; Persistency: " << energyCut
84  << " MeV; TofLimit: " << theTofLimit << " ns"
85  << "\n FrameRotation type " << rotType << " rTracker(cm)= " << rTracker / CLHEP::cm
86  << " zTracker(cm)= " << zTracker / CLHEP::cm
87  << " allowZeroEnergyLoss: " << allowZeroEnergyLoss
88  << " neverAccumulate: " << neverAccumulate << " printHits: " << printHits;
89 
90  slaveLowTof = std::make_unique<TrackingSlaveSD>(name + "LowTof");
91  slaveHighTof = std::make_unique<TrackingSlaveSD>(name + "HighTof");
92 
93  std::vector<std::string> temp;
94  temp.push_back(slaveLowTof.get()->name());
95  temp.push_back(slaveHighTof.get()->name());
96  setNames(temp);
97 
98  theG4ProcTypeEnumerator = std::make_unique<G4ProcessTypeEnumerator>();
99  theNumberingScheme = nullptr;
100 }
Log< level::Info, true > LogVerbatim
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::unique_ptr< TrackingSlaveSD > slaveHighTof
std::unique_ptr< FrameRotation > theRotation
std::unique_ptr< const G4ProcessTypeEnumerator > theG4ProcTypeEnumerator
std::unique_ptr< TrackingSlaveSD > slaveLowTof
double f[11][100]
void setNames(const std::vector< std::string > &)
SensitiveTkDetector(const std::string &iname, const SensitiveDetectorCatalog &clg)
TrackerG4SimHitNumberingScheme * theNumberingScheme

◆ ~TkAccumulatingSensitiveDetector()

TkAccumulatingSensitiveDetector::~TkAccumulatingSensitiveDetector ( )
override

Definition at line 102 of file TkAccumulatingSensitiveDetector.cc.

102 {}

Member Function Documentation

◆ clearHits()

void TkAccumulatingSensitiveDetector::clearHits ( )
overridevirtual

Implements SensitiveDetector.

Definition at line 354 of file TkAccumulatingSensitiveDetector.cc.

References slaveHighTof, and slaveLowTof.

Referenced by update().

354  {
355  slaveLowTof.get()->Initialize();
356  slaveHighTof.get()->Initialize();
357 }
std::unique_ptr< TrackingSlaveSD > slaveHighTof
std::unique_ptr< TrackingSlaveSD > slaveLowTof

◆ closeHit()

bool TkAccumulatingSensitiveDetector::closeHit ( const G4Step *  aStep)
private

Definition at line 331 of file TkAccumulatingSensitiveDetector.cc.

References PSimHit::exitPoint(), SensitiveDetector::LocalPreStepPosition(), LogDebug, mag(), mag2(), mySimHit, and theRotation.

Referenced by newHit().

331  {
332  const float tolerance2 = 0.0025f; // (0.5 mm)^2 are allowed between entry and exit
333  Local3DPoint theEntryPoint = theRotation.get()->transformPoint(LocalPreStepPosition(aStep));
334  LogDebug("TrackerSimDebug") << " closeHit: distance = " << (mySimHit->exitPoint() - theEntryPoint).mag();
335 
336  return ((mySimHit->exitPoint() - theEntryPoint).mag2() < tolerance2) ? true : false;
337 }
std::unique_ptr< FrameRotation > theRotation
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:46
Local3DPoint LocalPreStepPosition(const G4Step *step) const
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
T mag2() const
The vector magnitude squared. Equivalent to vec.dot(vec)
#define LogDebug(id)

◆ createHit()

void TkAccumulatingSensitiveDetector::createHit ( const G4Step *  aStep)
private

Definition at line 204 of file TkAccumulatingSensitiveDetector.cc.

References SensitiveDetector::cmsTrackInformation(), SensitiveDetector::ConvertToLocal3DPoint(), PSimHit::detUnitId(), PSimHit::energyLoss(), PSimHit::entryPoint(), Exception, PSimHit::exitPoint(), globalEntryPoint, globalExitPoint, lastId, lastTrack, SensitiveDetector::LocalPostStepPosition(), SensitiveDetector::LocalPreStepPosition(), LogDebug, mySimHit, oldVolume, G4TrackToParticleID::particleID(), PV3DBase< T, PVType, FrameType >::phi(), pname, printHits, px, py, pz, setDetUnitId(), groupFilesInBlocks::temp, theG4ProcTypeEnumerator, theRotation, PV3DBase< T, PVType, FrameType >::theta(), and PSimHit::trackId().

Referenced by ProcessHits().

204  {
205  // VI: previous hit should be already deleted
206  // in past here was a check if a hit is inside a sensitive detector,
207  // this is not needed, because call to senstive detector happens
208  // only inside the volume
209  const G4Track* theTrack = aStep->GetTrack();
210  Local3DPoint theExitPoint = theRotation.get()->transformPoint(LocalPostStepPosition(aStep));
211  Local3DPoint theEntryPoint;
212  //
213  // Check particle type - for gamma and neutral hadrons energy deposition
214  // should be local (VI)
215  //
216  if (0.0 == theTrack->GetDefinition()->GetPDGCharge()) {
217  theEntryPoint = theExitPoint;
218  } else {
219  theEntryPoint = theRotation.get()->transformPoint(LocalPreStepPosition(aStep));
220  }
221 
222  //
223  // This allows to send he skipEvent if it is outside!
224  //
225  const G4StepPoint* preStepPoint = aStep->GetPreStepPoint();
226  float thePabs = preStepPoint->GetMomentum().mag() / GeV;
227  float theTof = preStepPoint->GetGlobalTime() / nanosecond;
228  float theEnergyLoss = aStep->GetTotalEnergyDeposit() / GeV;
229  int theParticleType = G4TrackToParticleID::particleID(theTrack);
230  uint32_t theDetUnitId = setDetUnitId(aStep);
231  int theTrackID = theTrack->GetTrackID();
232  if (theDetUnitId == 0) {
233  edm::LogWarning("TkAccumulatingSensitiveDetector::createHit") << " theDetUnitId is not valid for " << GetName();
234  throw cms::Exception("TkAccumulatingSensitiveDetector::createHit")
235  << "cannot get theDetUnitId for G4Track " << theTrackID;
236  }
237 
238  // To whom assign the Hit?
239  // First iteration: if the track is to be stored, use the current number;
240  // otherwise, get to the mother
241  unsigned int theTrackIDInsideTheSimHit = theTrackID;
242 
243  const TrackInformation* temp = cmsTrackInformation(theTrack);
244  if (!temp->storeTrack()) {
245  // Go to the mother!
246  theTrackIDInsideTheSimHit = theTrack->GetParentID();
247  LogDebug("TrackerSimDebug") << " TkAccumulatingSensitiveDetector::createHit(): setting the TrackID from "
248  << theTrackIDInsideTheSimHit << " to the mother one " << theTrackIDInsideTheSimHit
249  << " " << theEnergyLoss;
250  } else {
251  LogDebug("TrackerSimDebug") << " TkAccumulatingSensitiveDetector:createHit(): leaving the current TrackID "
252  << theTrackIDInsideTheSimHit;
253  }
254 
255  const G4ThreeVector& gmd = preStepPoint->GetMomentumDirection();
256  // convert it to local frame
257  G4ThreeVector lmd =
258  ((G4TouchableHistory*)(preStepPoint->GetTouchable()))->GetHistory()->GetTopTransform().TransformAxis(gmd);
259  Local3DPoint lnmd = theRotation.get()->transformPoint(ConvertToLocal3DPoint(lmd));
260  float theThetaAtEntry = lnmd.theta();
261  float thePhiAtEntry = lnmd.phi();
262 
263  mySimHit = new UpdatablePSimHit(theEntryPoint,
264  theExitPoint,
265  thePabs,
266  theTof,
267  theEnergyLoss,
268  theParticleType,
269  theDetUnitId,
270  theTrackIDInsideTheSimHit,
271  theThetaAtEntry,
272  thePhiAtEntry,
273  theG4ProcTypeEnumerator.get()->processId(theTrack->GetCreatorProcess()));
274  lastId = theDetUnitId;
275  lastTrack = theTrackID;
276 
277  // only for debugging
278  if (printHits) {
279  // point on Geant4 unit (mm)
280  globalEntryPoint = ConvertToLocal3DPoint(preStepPoint->GetPosition());
281  globalExitPoint = ConvertToLocal3DPoint(aStep->GetPostStepPoint()->GetPosition());
282  // in CMS unit (GeV)
283  px = preStepPoint->GetMomentum().x() / CLHEP::GeV;
284  py = preStepPoint->GetMomentum().y() / CLHEP::GeV;
285  pz = preStepPoint->GetMomentum().z() / CLHEP::GeV;
286  oldVolume = preStepPoint->GetPhysicalVolume();
287  pname = theTrack->GetDefinition()->GetParticleName();
288  LogDebug("TrackerSimDebug") << " Created PSimHit: " << pname << " " << mySimHit->detUnitId() << " "
289  << mySimHit->trackId() << " " << theTrackID
290  << " p= " << aStep->GetPreStepPoint()->GetMomentum().mag() << " "
291  << mySimHit->energyLoss() << " " << mySimHit->entryPoint() << " "
292  << mySimHit->exitPoint();
293  }
294 }
std::unique_ptr< FrameRotation > theRotation
unsigned int detUnitId() const
Definition: PSimHit.h:97
std::unique_ptr< const G4ProcessTypeEnumerator > theG4ProcTypeEnumerator
Geom::Phi< T > phi() const
Definition: PV3DBase.h:66
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:46
uint32_t setDetUnitId(const G4Step *) override
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
Local3DPoint entryPoint() const
Entry point in the local Det frame.
Definition: PSimHit.h:43
unsigned int trackId() const
Definition: PSimHit.h:106
TrackInformation * cmsTrackInformation(const G4Track *aTrack)
static int particleID(const G4Track *)
Local3DPoint LocalPreStepPosition(const G4Step *step) const
float energyLoss() const
The energy deposit in the PSimHit, in ???.
Definition: PSimHit.h:79
Local3DPoint LocalPostStepPosition(const G4Step *step) const
Log< level::Warning, false > LogWarning
Geom::Theta< T > theta() const
Definition: PV3DBase.h:72
#define LogDebug(id)

◆ EndOfEvent()

void TkAccumulatingSensitiveDetector::EndOfEvent ( G4HCofThisEvent *  )
override

Definition at line 339 of file TkAccumulatingSensitiveDetector.cc.

References LogDebug, mySimHit, and sendHit().

339  {
340  LogDebug("TrackerSimDebug") << " Saving the last hit in a ROU " << GetName();
341  if (mySimHit != nullptr)
342  sendHit();
343 }
#define LogDebug(id)

◆ fillHits()

void TkAccumulatingSensitiveDetector::fillHits ( edm::PSimHitContainer cc,
const std::string &  hname 
)
overridevirtual

Implements SensitiveTkDetector.

Definition at line 359 of file TkAccumulatingSensitiveDetector.cc.

References gpuPixelDoublets::cc, slaveHighTof, and slaveLowTof.

359  {
360  if (slaveLowTof.get()->name() == hname) {
361  cc = slaveLowTof.get()->hits();
362  } else if (slaveHighTof.get()->name() == hname) {
363  cc = slaveHighTof.get()->hits();
364  }
365 }
std::unique_ptr< TrackingSlaveSD > slaveHighTof
uint32_t cc[maxCellsPerHit]
Definition: gpuFishbone.h:49
std::unique_ptr< TrackingSlaveSD > slaveLowTof

◆ newHit()

bool TkAccumulatingSensitiveDetector::newHit ( const G4Step *  aStep)
private

Definition at line 314 of file TkAccumulatingSensitiveDetector.cc.

References closeHit(), funct::false, lastId, lastTrack, LogDebug, and setDetUnitId().

Referenced by ProcessHits().

314  {
315  const G4Track* theTrack = aStep->GetTrack();
316 
317  // for neutral particles do not merge hits (V.I.)
318  if (0.0 == theTrack->GetDefinition()->GetPDGCharge())
319  return true;
320 
321  uint32_t theDetUnitId = setDetUnitId(aStep);
322  int theTrackID = theTrack->GetTrackID();
323 
324  LogDebug("TrackerSimDebug") << "newHit: OLD(detID,trID) = (" << lastId << "," << lastTrack << "), NEW = ("
325  << theDetUnitId << "," << theTrackID << ") Step length(mm)= " << aStep->GetStepLength()
326  << " Edep= " << aStep->GetTotalEnergyDeposit()
327  << " p= " << aStep->GetPreStepPoint()->GetMomentum().mag();
328  return ((theTrackID == lastTrack) && (lastId == theDetUnitId) && closeHit(aStep)) ? false : true;
329 }
uint32_t setDetUnitId(const G4Step *) override
#define LogDebug(id)

◆ ProcessHits()

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

Implements SensitiveDetector.

Definition at line 104 of file TkAccumulatingSensitiveDetector.cc.

References allowZeroEnergyLoss, createHit(), LogDebug, mySimHit, neverAccumulate, newHit(), sendHit(), and updateHit().

Referenced by LaserAlignmentSimulation::update().

104  {
105  LogDebug("TrackerSimDebug") << " Entering a new Step " << aStep->GetTotalEnergyDeposit() << " "
106  << aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName();
107 
108  if (aStep->GetTotalEnergyDeposit() > 0. || allowZeroEnergyLoss) {
109  if (!mySimHit) {
110  createHit(aStep);
111  } else if (neverAccumulate || newHit(aStep)) {
112  sendHit();
113  createHit(aStep);
114  } else {
115  updateHit(aStep);
116  }
117  return true;
118  }
119  return false;
120 }
#define LogDebug(id)

◆ sendHit()

void TkAccumulatingSensitiveDetector::sendHit ( )
private

Definition at line 171 of file TkAccumulatingSensitiveDetector.cc.

References PSimHit::detUnitId(), PSimHit::energyLoss(), PSimHit::entryPoint(), eventno, PSimHit::exitPoint(), globalEntryPoint, globalExitPoint, lastId, lastTrack, LogDebug, mySimHit, oldVolume, PSimHit::pabs(), pname, TkSimHitPrinter::printGlobal(), TkSimHitPrinter::printGlobalMomentum(), TkSimHitPrinter::printHitData(), printHits, TkSimHitPrinter::printLocal(), px, py, pz, slaveHighTof, slaveLowTof, TkSimHitPrinter::startNewSimHit(), theTofLimit, PSimHit::timeOfFlight(), and PSimHit::trackId().

Referenced by EndOfEvent(), and ProcessHits().

171  {
172  if (mySimHit == nullptr)
173  return;
174  if (printHits) {
175  TkSimHitPrinter thePrinter("TkHitPositionOSCAR.dat");
176  thePrinter.startNewSimHit(GetName(),
177  oldVolume->GetLogicalVolume()->GetName(),
178  mySimHit->detUnitId(),
179  mySimHit->trackId(),
180  lastTrack,
181  eventno);
182  thePrinter.printLocal(mySimHit->entryPoint(), mySimHit->exitPoint());
183  thePrinter.printGlobal(globalEntryPoint, globalExitPoint);
184  thePrinter.printHitData(pname, mySimHit->pabs(), mySimHit->energyLoss(), mySimHit->timeOfFlight());
185  thePrinter.printGlobalMomentum(px, py, pz);
186  LogDebug("TrackerSimDebug") << " Storing PSimHit: " << mySimHit->detUnitId() << " " << mySimHit->trackId() << " "
187  << mySimHit->energyLoss() << " " << mySimHit->entryPoint() << " "
188  << mySimHit->exitPoint();
189  }
190 
191  if (mySimHit->timeOfFlight() < theTofLimit) {
192  slaveLowTof.get()->processHits(*mySimHit); // implicit conversion (slicing) to PSimHit!!!
193  } else {
194  slaveHighTof.get()->processHits(*mySimHit); // implicit conversion (slicing) to PSimHit!!!
195  }
196  //
197  // clean up
198  delete mySimHit;
199  mySimHit = nullptr;
200  lastTrack = 0;
201  lastId = 0;
202 }
std::unique_ptr< TrackingSlaveSD > slaveHighTof
unsigned int detUnitId() const
Definition: PSimHit.h:97
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:46
std::unique_ptr< TrackingSlaveSD > slaveLowTof
Local3DPoint entryPoint() const
Entry point in the local Det frame.
Definition: PSimHit.h:43
unsigned int trackId() const
Definition: PSimHit.h:106
float pabs() const
fast and more accurate access to momentumAtEntry().mag()
Definition: PSimHit.h:67
float energyLoss() const
The energy deposit in the PSimHit, in ???.
Definition: PSimHit.h:79
float timeOfFlight() const
Definition: PSimHit.h:73
#define LogDebug(id)

◆ setDetUnitId()

uint32_t TkAccumulatingSensitiveDetector::setDetUnitId ( const G4Step *  step)
overridevirtual

Implements SensitiveDetector.

Definition at line 122 of file TkAccumulatingSensitiveDetector.cc.

References TrackerG4SimHitNumberingScheme::g4ToNumberingScheme(), and theNumberingScheme.

Referenced by createHit(), and newHit().

122  {
123  return theNumberingScheme->g4ToNumberingScheme(step->GetPreStepPoint()->GetTouchable());
124 }
unsigned int g4ToNumberingScheme(const G4VTouchable *)
step
Definition: StallMonitor.cc:98
TrackerG4SimHitNumberingScheme * theNumberingScheme

◆ update() [1/3]

void TkAccumulatingSensitiveDetector::update ( const BeginOfEvent )
overrideprotectedvirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfEvent *>.

Definition at line 345 of file TkAccumulatingSensitiveDetector.cc.

References clearHits(), eventno, and mySimHit.

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

345  {
346  clearHits();
347  eventno = (*i)()->GetEventID();
348  delete mySimHit;
349  mySimHit = nullptr;
350 }

◆ update() [2/3]

void TkAccumulatingSensitiveDetector::update ( const BeginOfTrack )
overrideprotectedvirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfTrack *>.

Definition at line 126 of file TkAccumulatingSensitiveDetector.cc.

References funct::abs(), SensitiveDetector::cmsTrackInformation(), energyCut, energyHistoryCut, info(), LogDebug, rTracker2, and zTracker.

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

126  {
127  const G4Track* gTrack = (*bot)();
128 
129 #ifdef DUMPPROCESSES
130  if (gTrack->GetCreatorProcess()) {
131  edm::LogVerbatim("TrackerSim") << " -> PROCESS CREATOR : " << gTrack->GetCreatorProcess()->GetProcessName();
132  } else {
133  edm::LogVerbatim("TrackerSim") << " -> No Creator process";
134  }
135 #endif
136 
137  //
138  //Position
139  //
140  const G4ThreeVector& pos = gTrack->GetPosition();
141  LogDebug("TrackerSimDebug") << " update(..) of " << gTrack->GetDefinition()->GetParticleName()
142  << " trackID= " << gTrack->GetTrackID() << " E(MeV)= " << gTrack->GetKineticEnergy()
143  << " Ecut= " << energyCut << " R(mm)= " << pos.perp() << " Z(mm)= " << pos.z();
144 
145  //
146  // Check if in Tracker Volume
147  //
148  if (pos.x() * pos.x() + pos.y() * pos.y() < rTracker2 && std::abs(pos.z()) < zTracker) {
149  //
150  // inside the Tracker
151  //
152  TrackInformation* info = nullptr;
153  if (gTrack->GetKineticEnergy() > energyCut) {
154  info = cmsTrackInformation(gTrack);
155  info->setStoreTrack();
156  }
157  //
158  // Save History?
159  //
160  if (gTrack->GetKineticEnergy() > energyHistoryCut) {
161  if (nullptr == info) {
162  info = cmsTrackInformation(gTrack);
163  }
164  info->putInHistory();
165  LogDebug("TrackerSimDebug") << " Track inside the tracker selected for HISTORY"
166  << " Track ID= " << gTrack->GetTrackID();
167  }
168  }
169 }
Log< level::Info, true > LogVerbatim
static const TGPicture * info(bool iBackgroundIsBlack)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
TrackInformation * cmsTrackInformation(const G4Track *aTrack)
#define LogDebug(id)

◆ update() [3/3]

void TkAccumulatingSensitiveDetector::update ( const BeginOfJob )
overrideprotectedvirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfJob *>.

Definition at line 352 of file TkAccumulatingSensitiveDetector.cc.

References numberingScheme(), pDD_, and theNumberingScheme.

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

static TrackerG4SimHitNumberingScheme & numberingScheme(const GeometricDet &det)
TrackerG4SimHitNumberingScheme * theNumberingScheme

◆ updateHit()

void TkAccumulatingSensitiveDetector::updateHit ( const G4Step *  aStep)
private

Definition at line 296 of file TkAccumulatingSensitiveDetector.cc.

References UpdatablePSimHit::addEnergyLoss(), SensitiveDetector::ConvertToLocal3DPoint(), PSimHit::detUnitId(), PSimHit::energyLoss(), PSimHit::entryPoint(), PSimHit::exitPoint(), globalExitPoint, SensitiveDetector::LocalPostStepPosition(), LogDebug, mySimHit, printHits, UpdatablePSimHit::setExitPoint(), theRotation, and PSimHit::trackId().

Referenced by ProcessHits().

296  {
297  // VI: in past here was a check if a hit is inside a sensitive detector,
298  // this is not needed, because call to senstive detector happens
299  // only inside the volume
300  Local3DPoint theExitPoint = theRotation.get()->transformPoint(LocalPostStepPosition(aStep));
301  float theEnergyLoss = aStep->GetTotalEnergyDeposit() / GeV;
302  mySimHit->setExitPoint(theExitPoint);
303  mySimHit->addEnergyLoss(theEnergyLoss);
304  if (printHits) {
305  globalExitPoint = ConvertToLocal3DPoint(aStep->GetPostStepPoint()->GetPosition());
306  LogDebug("TrackerSimDebug") << " updateHit: for " << aStep->GetTrack()->GetDefinition()->GetParticleName()
307  << " trackID= " << aStep->GetTrack()->GetTrackID() << " deltaEloss= " << theEnergyLoss
308  << "\n Updated PSimHit: " << mySimHit->detUnitId() << " " << mySimHit->trackId() << " "
309  << mySimHit->energyLoss() << " " << mySimHit->entryPoint() << " "
310  << mySimHit->exitPoint();
311  }
312 }
std::unique_ptr< FrameRotation > theRotation
unsigned int detUnitId() const
Definition: PSimHit.h:97
void setExitPoint(const Local3DPoint &exit)
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:46
Local3DPoint ConvertToLocal3DPoint(const G4ThreeVector &point) const
Local3DPoint entryPoint() const
Entry point in the local Det frame.
Definition: PSimHit.h:43
unsigned int trackId() const
Definition: PSimHit.h:106
void addEnergyLoss(float eloss)
float energyLoss() const
The energy deposit in the PSimHit, in ???.
Definition: PSimHit.h:79
Local3DPoint LocalPostStepPosition(const G4Step *step) const
#define LogDebug(id)

Member Data Documentation

◆ allowZeroEnergyLoss

bool TkAccumulatingSensitiveDetector::allowZeroEnergyLoss
private

◆ energyCut

float TkAccumulatingSensitiveDetector::energyCut
private

Definition at line 77 of file TkAccumulatingSensitiveDetector.h.

Referenced by TkAccumulatingSensitiveDetector(), and update().

◆ energyHistoryCut

float TkAccumulatingSensitiveDetector::energyHistoryCut
private

Definition at line 78 of file TkAccumulatingSensitiveDetector.h.

Referenced by TkAccumulatingSensitiveDetector(), and update().

◆ eventno

int TkAccumulatingSensitiveDetector::eventno
private

Definition at line 90 of file TkAccumulatingSensitiveDetector.h.

Referenced by sendHit(), and update().

◆ globalEntryPoint

Local3DPoint TkAccumulatingSensitiveDetector::globalEntryPoint
private

Definition at line 86 of file TkAccumulatingSensitiveDetector.h.

Referenced by createHit(), and sendHit().

◆ globalExitPoint

Local3DPoint TkAccumulatingSensitiveDetector::globalExitPoint
private

Definition at line 87 of file TkAccumulatingSensitiveDetector.h.

Referenced by createHit(), sendHit(), and updateHit().

◆ lastId

uint32_t TkAccumulatingSensitiveDetector::lastId
private

Definition at line 82 of file TkAccumulatingSensitiveDetector.h.

Referenced by createHit(), newHit(), and sendHit().

◆ lastTrack

int TkAccumulatingSensitiveDetector::lastTrack
private

Definition at line 83 of file TkAccumulatingSensitiveDetector.h.

Referenced by createHit(), newHit(), and sendHit().

◆ mySimHit

UpdatablePSimHit* TkAccumulatingSensitiveDetector::mySimHit
private

◆ neverAccumulate

bool TkAccumulatingSensitiveDetector::neverAccumulate
private

◆ oldVolume

const G4VPhysicalVolume* TkAccumulatingSensitiveDetector::oldVolume
private

Definition at line 88 of file TkAccumulatingSensitiveDetector.h.

Referenced by createHit(), and sendHit().

◆ pDD_

const GeometricDet* TkAccumulatingSensitiveDetector::pDD_
private

Definition at line 63 of file TkAccumulatingSensitiveDetector.h.

Referenced by update().

◆ pname

std::string TkAccumulatingSensitiveDetector::pname
private

Definition at line 91 of file TkAccumulatingSensitiveDetector.h.

Referenced by createHit(), and sendHit().

◆ printHits

bool TkAccumulatingSensitiveDetector::printHits
private

◆ px

float TkAccumulatingSensitiveDetector::px
private

Definition at line 89 of file TkAccumulatingSensitiveDetector.h.

Referenced by createHit(), and sendHit().

◆ py

float TkAccumulatingSensitiveDetector::py
private

Definition at line 89 of file TkAccumulatingSensitiveDetector.h.

Referenced by createHit(), and sendHit().

◆ pz

float TkAccumulatingSensitiveDetector::pz
private

Definition at line 89 of file TkAccumulatingSensitiveDetector.h.

Referenced by createHit(), and sendHit().

◆ rTracker

double TkAccumulatingSensitiveDetector::rTracker
private

Definition at line 74 of file TkAccumulatingSensitiveDetector.h.

Referenced by TkAccumulatingSensitiveDetector().

◆ rTracker2

double TkAccumulatingSensitiveDetector::rTracker2
private

Definition at line 73 of file TkAccumulatingSensitiveDetector.h.

Referenced by TkAccumulatingSensitiveDetector(), and update().

◆ slaveHighTof

std::unique_ptr<TrackingSlaveSD> TkAccumulatingSensitiveDetector::slaveHighTof
private

◆ slaveLowTof

std::unique_ptr<TrackingSlaveSD> TkAccumulatingSensitiveDetector::slaveLowTof
private

◆ theG4ProcTypeEnumerator

std::unique_ptr<const G4ProcessTypeEnumerator> TkAccumulatingSensitiveDetector::theG4ProcTypeEnumerator
private

Definition at line 68 of file TkAccumulatingSensitiveDetector.h.

Referenced by createHit(), and TkAccumulatingSensitiveDetector().

◆ theManager

const SimTrackManager* TkAccumulatingSensitiveDetector::theManager
private

Definition at line 64 of file TkAccumulatingSensitiveDetector.h.

◆ theNumberingScheme

TrackerG4SimHitNumberingScheme* TkAccumulatingSensitiveDetector::theNumberingScheme
private

◆ theRotation

std::unique_ptr<FrameRotation> TkAccumulatingSensitiveDetector::theRotation
private

◆ theTofLimit

float TkAccumulatingSensitiveDetector::theTofLimit
private

Definition at line 76 of file TkAccumulatingSensitiveDetector.h.

Referenced by sendHit(), and TkAccumulatingSensitiveDetector().

◆ zTracker

double TkAccumulatingSensitiveDetector::zTracker
private

Definition at line 75 of file TkAccumulatingSensitiveDetector.h.

Referenced by TkAccumulatingSensitiveDetector(), and update().