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 DDCompactView &, const SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
 
 ~TkAccumulatingSensitiveDetector () override
 
- Public Member Functions inherited from SensitiveTkDetector
 SensitiveTkDetector (const std::string &iname, const DDCompactView &cpv, const SensitiveDetectorCatalog &clg, edm::ParameterSet const &p)
 
- Public Member Functions inherited from SensitiveDetector
void EndOfEvent (G4HCofThisEvent *eventHC) override
 
const std::vector< std::string > & getNames () const
 
void Initialize (G4HCofThisEvent *eventHC) override
 
 SensitiveDetector (const std::string &iname, const DDCompactView &cpv, const SensitiveDetectorCatalog &, edm::ParameterSet const &p)
 
 ~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
 
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
 
std::unique_ptr< 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 28 of file TkAccumulatingSensitiveDetector.h.

Constructor & Destructor Documentation

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

Definition at line 48 of file TkAccumulatingSensitiveDetector.cc.

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

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

Definition at line 99 of file TkAccumulatingSensitiveDetector.cc.

100 {}

Member Function Documentation

void TkAccumulatingSensitiveDetector::clearHits ( )
overridevirtual

Implements SensitiveDetector.

Definition at line 375 of file TkAccumulatingSensitiveDetector.cc.

References slaveHighTof, and slaveLowTof.

Referenced by update().

376 {
377  slaveLowTof.get()->Initialize();
378  slaveHighTof.get()->Initialize();
379 }
std::unique_ptr< TrackingSlaveSD > slaveHighTof
std::unique_ptr< TrackingSlaveSD > slaveLowTof
bool TkAccumulatingSensitiveDetector::closeHit ( const G4Step *  aStep)
private

Definition at line 339 of file TkAccumulatingSensitiveDetector.cc.

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

Referenced by newHit().

340 {
341  const float tolerance2 = 0.0025f; // (0.5 mm)^2 are allowed between entry and exit
342  Local3DPoint theEntryPoint = theRotation.get()->transformPoint(LocalPreStepPosition(aStep));
343  LogDebug("TrackerSimDebug")<< " closeHit: distance = "
344  << (mySimHit->exitPoint()-theEntryPoint).mag();
345 
346  return ((mySimHit->exitPoint()-theEntryPoint).mag2() < tolerance2) ? true : false;
347 }
#define LogDebug(id)
std::unique_ptr< FrameRotation > theRotation
T mag() const
The vector magnitude. Equivalent to sqrt(vec.mag2())
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:38
T mag2() const
The vector magnitude squared. Equivalent to vec.dot(vec)
Local3DPoint LocalPreStepPosition(const G4Step *step) const
void TkAccumulatingSensitiveDetector::createHit ( const G4Step *  aStep)
private

Definition at line 209 of file TkAccumulatingSensitiveDetector.cc.

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

Referenced by ProcessHits().

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

Definition at line 349 of file TkAccumulatingSensitiveDetector.cc.

References LogDebug, mySimHit, and sendHit().

350 {
351  LogDebug("TrackerSimDebug")<< " Saving the last hit in a ROU " << GetName();
352  if (mySimHit != nullptr) sendHit();
353 }
#define LogDebug(id)
void TkAccumulatingSensitiveDetector::fillHits ( edm::PSimHitContainer cc,
const std::string &  hname 
)
overridevirtual

Implements SensitiveTkDetector.

Definition at line 381 of file TkAccumulatingSensitiveDetector.cc.

References slaveHighTof, and slaveLowTof.

381  {
382 
383  if (slaveLowTof.get()->name() == hname) { cc=slaveLowTof.get()->hits(); }
384  else if (slaveHighTof.get()->name() == hname) { cc=slaveHighTof.get()->hits();}
385 }
std::unique_ptr< TrackingSlaveSD > slaveHighTof
std::unique_ptr< TrackingSlaveSD > slaveLowTof
bool TkAccumulatingSensitiveDetector::newHit ( const G4Step *  aStep)
private

Definition at line 320 of file TkAccumulatingSensitiveDetector.cc.

References closeHit(), lastId, lastTrack, LogDebug, and setDetUnitId().

Referenced by ProcessHits().

321 {
322  const G4Track * theTrack = aStep->GetTrack();
323 
324  // for neutral particles do not merge hits (V.I.)
325  if(0.0 == theTrack->GetDefinition()->GetPDGCharge()) return true;
326 
327  uint32_t theDetUnitId = setDetUnitId(aStep);
328  int theTrackID = theTrack->GetTrackID();
329 
330  LogDebug("TrackerSimDebug")<< "newHit: OLD(detID,trID) = (" << lastId << "," << lastTrack
331  << "), NEW = (" << theDetUnitId << "," << theTrackID
332  << ") Step length(mm)= " << aStep->GetStepLength()
333  << " Edep= " << aStep->GetTotalEnergyDeposit()
334  << " p= " << aStep->GetPreStepPoint()->GetMomentum().mag();
335  return ((theTrackID == lastTrack) && (lastId == theDetUnitId) && closeHit(aStep))
336  ? false : true;
337 }
#define LogDebug(id)
uint32_t setDetUnitId(const G4Step *) override
bool TkAccumulatingSensitiveDetector::ProcessHits ( G4Step *  aStep,
G4TouchableHistory *   
)
overridevirtual

Implements SensitiveDetector.

Definition at line 102 of file TkAccumulatingSensitiveDetector.cc.

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

Referenced by LaserAlignmentSimulation::update().

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

Definition at line 179 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().

180 {
181  if (mySimHit == nullptr) return;
182  if (printHits)
183  {
184  TkSimHitPrinter thePrinter("TkHitPositionOSCAR.dat");
185  thePrinter.startNewSimHit(GetName(),oldVolume->GetLogicalVolume()->GetName(),
187  thePrinter.printLocal(mySimHit->entryPoint(),mySimHit->exitPoint());
188  thePrinter.printGlobal(globalEntryPoint,globalExitPoint);
189  thePrinter.printHitData(pname,mySimHit->pabs(),mySimHit->energyLoss(),mySimHit->timeOfFlight());
190  thePrinter.printGlobalMomentum(px,py,pz);
191  LogDebug("TrackerSimDebug")<< " Storing PSimHit: " << mySimHit->detUnitId()
192  << " " << mySimHit->trackId() << " " << mySimHit->energyLoss()
193  << " " << mySimHit->entryPoint() << " " << mySimHit->exitPoint();
194  }
195 
196  if (mySimHit->timeOfFlight() < theTofLimit) {
197  slaveLowTof.get()->processHits(*mySimHit); // implicit conversion (slicing) to PSimHit!!!
198  } else {
199  slaveHighTof.get()->processHits(*mySimHit); // implicit conversion (slicing) to PSimHit!!!
200  }
201  //
202  // clean up
203  delete mySimHit;
204  mySimHit = nullptr;
205  lastTrack = 0;
206  lastId = 0;
207 }
#define LogDebug(id)
std::unique_ptr< TrackingSlaveSD > slaveHighTof
std::unique_ptr< TrackingSlaveSD > slaveLowTof
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:38
float timeOfFlight() const
Definition: PSimHit.h:69
float pabs() const
fast and more accurate access to momentumAtEntry().mag()
Definition: PSimHit.h:63
float energyLoss() const
The energy deposit in the PSimHit, in ???.
Definition: PSimHit.h:75
unsigned int trackId() const
Definition: PSimHit.h:102
Local3DPoint entryPoint() const
Entry point in the local Det frame.
Definition: PSimHit.h:35
unsigned int detUnitId() const
Definition: PSimHit.h:93
uint32_t TkAccumulatingSensitiveDetector::setDetUnitId ( const G4Step *  step)
overridevirtual

Implements SensitiveDetector.

Definition at line 127 of file TkAccumulatingSensitiveDetector.cc.

References theNumberingScheme.

Referenced by createHit(), and newHit().

128 {
129  return theNumberingScheme.get()->g4ToNumberingScheme(step->GetPreStepPoint()->GetTouchable());
130 }
std::unique_ptr< TrackerG4SimHitNumberingScheme > theNumberingScheme
step
void TkAccumulatingSensitiveDetector::update ( const BeginOfEvent )
overrideprotectedvirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfEvent * >.

Definition at line 355 of file TkAccumulatingSensitiveDetector.cc.

References clearHits(), eventno, and mySimHit.

Referenced by progressbar.ProgressBar::__next__(), MatrixUtil.Matrix::__setitem__(), MatrixUtil.Steps::__setitem__(), Vispa.Gui.VispaWidget.VispaWidget::autosize(), Vispa.Views.LineDecayView.LineDecayContainer::createObject(), Vispa.Views.LineDecayView.LineDecayContainer::deselectAllObjects(), Vispa.Gui.VispaWidgetOwner.VispaWidgetOwner::deselectAllWidgets(), Vispa.Gui.VispaWidget.VispaWidget::enableAutosizing(), progressbar.ProgressBar::finish(), Vispa.Gui.MenuWidget.MenuWidget::leaveEvent(), Vispa.Gui.VispaWidgetOwner.VispaWidgetOwner::mouseMoveEvent(), Vispa.Gui.MenuWidget.MenuWidget::mouseMoveEvent(), Vispa.Views.LineDecayView.LineDecayContainer::mouseMoveEvent(), Vispa.Gui.VispaWidgetOwner.VispaWidgetOwner::mouseReleaseEvent(), Vispa.Views.LineDecayView.LineDecayContainer::objectMoved(), MatrixUtil.Steps::overwrite(), Vispa.Views.LineDecayView.LineDecayContainer::removeObject(), Vispa.Gui.ConnectableWidget.ConnectableWidget::removePorts(), Vispa.Gui.FindDialog.FindDialog::reset(), Vispa.Gui.PortConnection.PointToPointConnection::select(), Vispa.Gui.VispaWidget.VispaWidget::select(), Vispa.Views.LineDecayView.LineDecayContainer::select(), Vispa.Gui.VispaWidget.VispaWidget::setText(), Vispa.Gui.VispaWidget.VispaWidget::setTitle(), Vispa.Gui.ZoomableWidget.ZoomableWidget::setZoom(), Vispa.Views.LineDecayView.LineDecayContainer::setZoom(), and Vispa.Gui.PortConnection.PointToPointConnection::updateConnection().

356 {
357  clearHits();
358  eventno = (*i)()->GetEventID();
359  delete mySimHit;
360  mySimHit = nullptr;
361 }
void TkAccumulatingSensitiveDetector::update ( const BeginOfTrack )
overrideprotectedvirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfTrack * >.

Definition at line 132 of file TkAccumulatingSensitiveDetector.cc.

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

Referenced by progressbar.ProgressBar::__next__(), MatrixUtil.Matrix::__setitem__(), MatrixUtil.Steps::__setitem__(), Vispa.Gui.VispaWidget.VispaWidget::autosize(), Vispa.Views.LineDecayView.LineDecayContainer::createObject(), Vispa.Views.LineDecayView.LineDecayContainer::deselectAllObjects(), Vispa.Gui.VispaWidgetOwner.VispaWidgetOwner::deselectAllWidgets(), Vispa.Gui.VispaWidget.VispaWidget::enableAutosizing(), progressbar.ProgressBar::finish(), Vispa.Gui.MenuWidget.MenuWidget::leaveEvent(), Vispa.Gui.VispaWidgetOwner.VispaWidgetOwner::mouseMoveEvent(), Vispa.Gui.MenuWidget.MenuWidget::mouseMoveEvent(), Vispa.Views.LineDecayView.LineDecayContainer::mouseMoveEvent(), Vispa.Gui.VispaWidgetOwner.VispaWidgetOwner::mouseReleaseEvent(), Vispa.Views.LineDecayView.LineDecayContainer::objectMoved(), MatrixUtil.Steps::overwrite(), Vispa.Views.LineDecayView.LineDecayContainer::removeObject(), Vispa.Gui.ConnectableWidget.ConnectableWidget::removePorts(), Vispa.Gui.FindDialog.FindDialog::reset(), Vispa.Gui.PortConnection.PointToPointConnection::select(), Vispa.Gui.VispaWidget.VispaWidget::select(), Vispa.Views.LineDecayView.LineDecayContainer::select(), Vispa.Gui.VispaWidget.VispaWidget::setText(), Vispa.Gui.VispaWidget.VispaWidget::setTitle(), Vispa.Gui.ZoomableWidget.ZoomableWidget::setZoom(), Vispa.Views.LineDecayView.LineDecayContainer::setZoom(), and Vispa.Gui.PortConnection.PointToPointConnection::updateConnection().

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

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfJob * >.

Definition at line 363 of file TkAccumulatingSensitiveDetector.cc.

References edm::EventSetup::get(), numberingScheme(), and theNumberingScheme.

Referenced by progressbar.ProgressBar::__next__(), MatrixUtil.Matrix::__setitem__(), MatrixUtil.Steps::__setitem__(), Vispa.Gui.VispaWidget.VispaWidget::autosize(), Vispa.Views.LineDecayView.LineDecayContainer::createObject(), Vispa.Views.LineDecayView.LineDecayContainer::deselectAllObjects(), Vispa.Gui.VispaWidgetOwner.VispaWidgetOwner::deselectAllWidgets(), Vispa.Gui.VispaWidget.VispaWidget::enableAutosizing(), progressbar.ProgressBar::finish(), Vispa.Gui.MenuWidget.MenuWidget::leaveEvent(), Vispa.Gui.VispaWidgetOwner.VispaWidgetOwner::mouseMoveEvent(), Vispa.Gui.MenuWidget.MenuWidget::mouseMoveEvent(), Vispa.Views.LineDecayView.LineDecayContainer::mouseMoveEvent(), Vispa.Gui.VispaWidgetOwner.VispaWidgetOwner::mouseReleaseEvent(), Vispa.Views.LineDecayView.LineDecayContainer::objectMoved(), MatrixUtil.Steps::overwrite(), Vispa.Views.LineDecayView.LineDecayContainer::removeObject(), Vispa.Gui.ConnectableWidget.ConnectableWidget::removePorts(), Vispa.Gui.FindDialog.FindDialog::reset(), Vispa.Gui.PortConnection.PointToPointConnection::select(), Vispa.Gui.VispaWidget.VispaWidget::select(), Vispa.Views.LineDecayView.LineDecayContainer::select(), Vispa.Gui.VispaWidget.VispaWidget::setText(), Vispa.Gui.VispaWidget.VispaWidget::setTitle(), Vispa.Gui.ZoomableWidget.ZoomableWidget::setZoom(), Vispa.Views.LineDecayView.LineDecayContainer::setZoom(), and Vispa.Gui.PortConnection.PointToPointConnection::updateConnection().

364 {
366  const edm::EventSetup* es = (*i)();
367  es->get<IdealGeometryRecord>().get( pDD );
368 
370  es->get<IdealGeometryRecord>().get(pView);
371 
372  theNumberingScheme.reset(&(numberingScheme(*pView,*pDD)));
373 }
static TrackerG4SimHitNumberingScheme & numberingScheme(const DDCompactView &cpv, const GeometricDet &det)
std::unique_ptr< TrackerG4SimHitNumberingScheme > theNumberingScheme
T get() const
Definition: EventSetup.h:68
void TkAccumulatingSensitiveDetector::updateHit ( const G4Step *  aStep)
private

Definition at line 300 of file TkAccumulatingSensitiveDetector.cc.

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

Referenced by ProcessHits().

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

Member Data Documentation

bool TkAccumulatingSensitiveDetector::allowZeroEnergyLoss
private
float TkAccumulatingSensitiveDetector::energyCut
private

Definition at line 75 of file TkAccumulatingSensitiveDetector.h.

Referenced by TkAccumulatingSensitiveDetector(), and update().

float TkAccumulatingSensitiveDetector::energyHistoryCut
private

Definition at line 76 of file TkAccumulatingSensitiveDetector.h.

Referenced by TkAccumulatingSensitiveDetector(), and update().

int TkAccumulatingSensitiveDetector::eventno
private

Definition at line 88 of file TkAccumulatingSensitiveDetector.h.

Referenced by sendHit(), and update().

Local3DPoint TkAccumulatingSensitiveDetector::globalEntryPoint
private

Definition at line 84 of file TkAccumulatingSensitiveDetector.h.

Referenced by createHit(), and sendHit().

Local3DPoint TkAccumulatingSensitiveDetector::globalExitPoint
private

Definition at line 85 of file TkAccumulatingSensitiveDetector.h.

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

uint32_t TkAccumulatingSensitiveDetector::lastId
private

Definition at line 80 of file TkAccumulatingSensitiveDetector.h.

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

int TkAccumulatingSensitiveDetector::lastTrack
private

Definition at line 81 of file TkAccumulatingSensitiveDetector.h.

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

UpdatablePSimHit* TkAccumulatingSensitiveDetector::mySimHit
private
bool TkAccumulatingSensitiveDetector::neverAccumulate
private
const G4VPhysicalVolume* TkAccumulatingSensitiveDetector::oldVolume
private

Definition at line 86 of file TkAccumulatingSensitiveDetector.h.

Referenced by createHit(), and sendHit().

std::string TkAccumulatingSensitiveDetector::pname
private

Definition at line 89 of file TkAccumulatingSensitiveDetector.h.

Referenced by createHit(), and sendHit().

bool TkAccumulatingSensitiveDetector::printHits
private
float TkAccumulatingSensitiveDetector::px
private

Definition at line 87 of file TkAccumulatingSensitiveDetector.h.

Referenced by createHit(), and sendHit().

float TkAccumulatingSensitiveDetector::py
private

Definition at line 87 of file TkAccumulatingSensitiveDetector.h.

Referenced by createHit(), and sendHit().

float TkAccumulatingSensitiveDetector::pz
private

Definition at line 87 of file TkAccumulatingSensitiveDetector.h.

Referenced by createHit(), and sendHit().

double TkAccumulatingSensitiveDetector::rTracker
private

Definition at line 72 of file TkAccumulatingSensitiveDetector.h.

Referenced by TkAccumulatingSensitiveDetector().

double TkAccumulatingSensitiveDetector::rTracker2
private

Definition at line 71 of file TkAccumulatingSensitiveDetector.h.

Referenced by TkAccumulatingSensitiveDetector(), and update().

std::unique_ptr<TrackingSlaveSD> TkAccumulatingSensitiveDetector::slaveHighTof
private
std::unique_ptr<TrackingSlaveSD> TkAccumulatingSensitiveDetector::slaveLowTof
private
std::unique_ptr<const G4ProcessTypeEnumerator> TkAccumulatingSensitiveDetector::theG4ProcTypeEnumerator
private

Definition at line 66 of file TkAccumulatingSensitiveDetector.h.

Referenced by createHit(), and TkAccumulatingSensitiveDetector().

const SimTrackManager* TkAccumulatingSensitiveDetector::theManager
private

Definition at line 62 of file TkAccumulatingSensitiveDetector.h.

std::unique_ptr<TrackerG4SimHitNumberingScheme> TkAccumulatingSensitiveDetector::theNumberingScheme
private
std::unique_ptr<FrameRotation> TkAccumulatingSensitiveDetector::theRotation
private
float TkAccumulatingSensitiveDetector::theTofLimit
private

Definition at line 74 of file TkAccumulatingSensitiveDetector.h.

Referenced by sendHit(), and TkAccumulatingSensitiveDetector().

double TkAccumulatingSensitiveDetector::zTracker
private

Definition at line 73 of file TkAccumulatingSensitiveDetector.h.

Referenced by TkAccumulatingSensitiveDetector(), and update().