CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
PLTSensitiveDetector Class Reference

#include <PLTSensitiveDetector.h>

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

Public Member Functions

virtual void EndOfEvent (G4HCofThisEvent *)
 
void fillHits (edm::PSimHitContainer &, std::string use)
 
 PLTSensitiveDetector (std::string, const DDCompactView &, SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *)
 
virtual bool ProcessHits (G4Step *, G4TouchableHistory *)
 
virtual uint32_t setDetUnitId (G4Step *)
 
virtual ~PLTSensitiveDetector ()
 
- Public Member Functions inherited from SensitiveTkDetector
 SensitiveTkDetector (std::string &iname, const DDCompactView &cpv, SensitiveDetectorCatalog &clg, edm::ParameterSet const &p)
 
- Public Member Functions inherited from SensitiveDetector
virtual void AssignSD (std::string &vname)
 
Local3DPoint ConvertToLocal3DPoint (G4ThreeVector point)
 
Local3DPoint FinalStepPosition (G4Step *s, coordinates)
 
virtual std::vector< std::string > getNames ()
 
virtual void Initialize (G4HCofThisEvent *eventHC)
 
Local3DPoint InitialStepPosition (G4Step *s, coordinates)
 
std::string nameOfSD ()
 
void NaNTrap (G4Step *step)
 
void Register ()
 
 SensitiveDetector (std::string &iname, const DDCompactView &cpv, SensitiveDetectorCatalog &, edm::ParameterSet const &p)
 
virtual ~SensitiveDetector ()
 
- 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 ()
 

Private Member Functions

virtual void clearHits ()
 
virtual bool closeHit (G4Step *)
 
virtual void createHit (G4Step *)
 
TrackInformationgetOrCreateTrackInformation (const G4Track *)
 
virtual bool newHit (G4Step *)
 
virtual void sendHit ()
 
void update (const BeginOfEvent *)
 This routine will be called when the appropriate signal arrives. More...
 
void update (const BeginOfTrack *)
 This routine will be called when the appropriate signal arrives. More...
 
void update (const BeginOfJob *)
 This routine will be called when the appropriate signal arrives. More...
 
virtual void updateHit (G4Step *)
 

Private Attributes

float energyCut
 
float energyHistoryCut
 
int eventno
 
Local3DPoint globalEntryPoint
 
Local3DPoint globalExitPoint
 
uint32_t lastId
 
unsigned int lastTrack
 
G4TrackToParticleIDmyG4TrackToParticleID
 
std::string myName
 
UpdatablePSimHitmySimHit
 
G4VPhysicalVolume * oldVolume
 
std::string pname
 
TrackingSlaveSDslave
 
G4ProcessTypeEnumeratortheG4ProcessTypeEnumerator
 

Additional Inherited Members

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

Detailed Description

Definition at line 30 of file PLTSensitiveDetector.h.

Constructor & Destructor Documentation

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

Definition at line 32 of file PLTSensitiveDetector.cc.

References SensitiveDetector::AssignSD(), energyCut, energyHistoryCut, edm::ParameterSet::getParameter(), SensitiveDetectorCatalog::logicalNames(), myG4TrackToParticleID, SensitiveDetector::Register(), slave, and theG4ProcessTypeEnumerator.

36  :
37  SensitiveTkDetector(name, cpv, clg, p), myName(name), mySimHit(0),
38  oldVolume(0), lastId(0), lastTrack(0), eventno(0) {
39 
40  edm::ParameterSet m_TrackerSD = p.getParameter<edm::ParameterSet>("PltSD");
41  energyCut = m_TrackerSD.getParameter<double>("EnergyThresholdForPersistencyInGeV")*GeV; //default must be 0.5
42  energyHistoryCut = m_TrackerSD.getParameter<double>("EnergyThresholdForHistoryInGeV")*GeV;//default must be 0.05
43 
44  edm::LogInfo("PLTSensitiveDetector") <<"Criteria for Saving Tracker SimTracks: \n "
45  <<" History: "<<energyHistoryCut<< " MeV ; Persistency: "<< energyCut<<" MeV\n"
46  <<" Constructing a PLTSensitiveDetector with ";
47 
48  slave = new TrackingSlaveSD(name);
49 
50  // Now attach the right detectors (LogicalVolumes) to me
51  std::vector<std::string> lvNames = clg.logicalNames(name);
52  this->Register();
53  for (std::vector<std::string>::iterator it = lvNames.begin();
54  it != lvNames.end(); it++) {
55  edm::LogInfo("PLTSensitiveDetector")<< name << " attaching LV " << *it;
56  this->AssignSD(*it);
57  }
58 
61 }
T getParameter(std::string const &) const
UpdatablePSimHit * mySimHit
std::vector< std::string > logicalNames(std::string &readoutName)
G4TrackToParticleID * myG4TrackToParticleID
G4ProcessTypeEnumerator * theG4ProcessTypeEnumerator
G4VPhysicalVolume * oldVolume
SensitiveTkDetector(std::string &iname, const DDCompactView &cpv, SensitiveDetectorCatalog &clg, edm::ParameterSet const &p)
TrackingSlaveSD * slave
virtual void AssignSD(std::string &vname)
PLTSensitiveDetector::~PLTSensitiveDetector ( )
virtual

Definition at line 63 of file PLTSensitiveDetector.cc.

References myG4TrackToParticleID, slave, and theG4ProcessTypeEnumerator.

63  {
64  delete slave;
66  delete myG4TrackToParticleID;
67 }
G4TrackToParticleID * myG4TrackToParticleID
G4ProcessTypeEnumerator * theG4ProcessTypeEnumerator
TrackingSlaveSD * slave

Member Function Documentation

void PLTSensitiveDetector::clearHits ( )
privatevirtual

Implements SensitiveDetector.

Definition at line 240 of file PLTSensitiveDetector.cc.

References TrackingSlaveSD::Initialize(), and slave.

Referenced by update().

240  {
241  slave->Initialize();
242 }
virtual void Initialize()
TrackingSlaveSD * slave
bool PLTSensitiveDetector::closeHit ( G4Step *  aStep)
privatevirtual

Definition at line 167 of file PLTSensitiveDetector.cc.

References PSimHit::exitPoint(), SensitiveDetector::InitialStepPosition(), SensitiveDetector::LocalCoordinates, LogDebug, mag(), and mySimHit.

Referenced by newHit().

167  {
168 
169  if (mySimHit == 0) return false;
170  const float tolerance = 0.05 * mm; // 50 micron are allowed between the exit
171  // point of the current hit and the entry point of the new hit
173  LogDebug("PLTSensitiveDetector")<< " closeHit: distance = " << (mySimHit->exitPoint()-theEntryPoint).mag();
174 
175  if ((mySimHit->exitPoint()-theEntryPoint).mag()<tolerance) return true;
176  return false;
177 }
#define LogDebug(id)
UpdatablePSimHit * mySimHit
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
Local3DPoint InitialStepPosition(G4Step *s, coordinates)
void PLTSensitiveDetector::createHit ( G4Step *  aStep)
privatevirtual

Definition at line 179 of file PLTSensitiveDetector.cc.

References SensitiveDetector::ConvertToLocal3DPoint(), PSimHit::detUnitId(), PSimHit::energyLoss(), PSimHit::entryPoint(), PSimHit::exitPoint(), SensitiveDetector::FinalStepPosition(), globalEntryPoint, globalExitPoint, SensitiveDetector::InitialStepPosition(), lastId, lastTrack, SensitiveDetector::LocalCoordinates, LogDebug, myG4TrackToParticleID, mySimHit, oldVolume, G4TrackToParticleID::particleID(), PV3DBase< T, PVType, FrameType >::phi(), pname, G4ProcessTypeEnumerator::processId(), setDetUnitId(), theG4ProcessTypeEnumerator, PV3DBase< T, PVType, FrameType >::theta(), PSimHit::trackId(), findQualityFiles::v, and SensitiveDetector::WorldCoordinates.

Referenced by ProcessHits().

179  {
180 
181  if (mySimHit != 0) {
182  delete mySimHit;
183  mySimHit=0;
184  }
185 
186  G4Track * theTrack = aStep->GetTrack();
187  G4VPhysicalVolume * v = aStep->GetPreStepPoint()->GetPhysicalVolume();
188 
191 
192  float thePabs = aStep->GetPreStepPoint()->GetMomentum().mag()/GeV;
193  float theTof = aStep->GetPreStepPoint()->GetGlobalTime()/nanosecond;
194  float theEnergyLoss = aStep->GetTotalEnergyDeposit()/GeV;
195  int theParticleType = myG4TrackToParticleID->particleID(theTrack);
196  uint32_t theDetUnitId = setDetUnitId(aStep);
197 
200  pname = theTrack->GetDynamicParticle()->GetDefinition()->GetParticleName();
201 
202  unsigned int theTrackID = theTrack->GetTrackID();
203 
204  G4ThreeVector gmd = aStep->GetPreStepPoint()->GetMomentumDirection();
205  // convert it to local frame
206  G4ThreeVector lmd = ((G4TouchableHistory *)(aStep->GetPreStepPoint()->GetTouchable()))->GetHistory()->GetTopTransform().TransformAxis(gmd);
208  float theThetaAtEntry = lnmd.theta();
209  float thePhiAtEntry = lnmd.phi();
210 
211  mySimHit = new UpdatablePSimHit(theEntryPoint,theExitPoint,thePabs,theTof,
212  theEnergyLoss,theParticleType,theDetUnitId,
213  theTrackID,theThetaAtEntry,thePhiAtEntry,
214  theG4ProcessTypeEnumerator->processId(theTrack->GetCreatorProcess()));
215 
216  LogDebug("PLTSensitiveDetector") << " Created PSimHit: " << pname << " "
217  << mySimHit->detUnitId() << " " << mySimHit->trackId()
218  << " " << mySimHit->energyLoss() << " "
219  << mySimHit->entryPoint() << " " << mySimHit->exitPoint();
220  lastId = theDetUnitId;
221  lastTrack = theTrackID;
222  oldVolume = v;
223 }
#define LogDebug(id)
UpdatablePSimHit * mySimHit
G4TrackToParticleID * myG4TrackToParticleID
Local3DPoint ConvertToLocal3DPoint(G4ThreeVector point)
Geom::Phi< T > phi() const
Definition: PV3DBase.h:69
G4ProcessTypeEnumerator * theG4ProcessTypeEnumerator
G4VPhysicalVolume * oldVolume
Geom::Theta< T > theta() const
Definition: PV3DBase.h:75
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:38
int particleID(const G4Track *)
virtual uint32_t setDetUnitId(G4Step *)
unsigned int processId(const G4VProcess *)
float energyLoss() const
The energy deposit in the PSimHit, in ???.
Definition: PSimHit.h:75
unsigned int trackId() const
Definition: PSimHit.h:102
Local3DPoint FinalStepPosition(G4Step *s, coordinates)
Local3DPoint entryPoint() const
Entry point in the local Det frame.
Definition: PSimHit.h:35
Local3DPoint InitialStepPosition(G4Step *s, coordinates)
unsigned int detUnitId() const
Definition: PSimHit.h:93
void PLTSensitiveDetector::EndOfEvent ( G4HCofThisEvent *  )
virtual

Reimplemented from SensitiveDetector.

Definition at line 109 of file PLTSensitiveDetector.cc.

References LogDebug, myName, mySimHit, and sendHit().

109  {
110 
111  LogDebug("PLTSensitiveDetector")<< " Saving the last hit in a ROU " << myName;
112 
113  if (mySimHit == 0) return;
114  sendHit();
115 }
#define LogDebug(id)
UpdatablePSimHit * mySimHit
void PLTSensitiveDetector::fillHits ( edm::PSimHitContainer c,
std::string  use 
)
virtual

Implements SensitiveTkDetector.

Definition at line 117 of file PLTSensitiveDetector.cc.

References TrackingSlaveSD::hits(), n, TrackingSlaveSD::name(), and slave.

117  {
118  if (slave->name() == n) c=slave->hits();
119 }
std::string name() const
std::vector< PSimHit > & hits()
TrackingSlaveSD * slave
TrackInformation * PLTSensitiveDetector::getOrCreateTrackInformation ( const G4Track *  gTrack)
private

Definition at line 244 of file PLTSensitiveDetector.cc.

References info, and groupFilesInBlocks::temp.

Referenced by ProcessHits().

244  {
245  G4VUserTrackInformation* temp = gTrack->GetUserInformation();
246  if (temp == 0){
247  edm::LogError("PLTSensitiveDetector") <<" ERROR: no G4VUserTrackInformation available";
248  abort();
249  }else{
250  TrackInformation* info = dynamic_cast<TrackInformation*>(temp);
251  if (info == 0){
252  edm::LogError("PLTSensitiveDetector") <<" ERROR: TkSimTrackSelection: the UserInformation does not appear to be a TrackInformation";
253  abort();
254  }
255  return info;
256  }
257 }
bool PLTSensitiveDetector::newHit ( G4Step *  aStep)
privatevirtual

Definition at line 153 of file PLTSensitiveDetector.cc.

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

Referenced by ProcessHits().

153  {
154 
155  G4Track * theTrack = aStep->GetTrack();
156  uint32_t theDetUnitId = setDetUnitId(aStep);
157  unsigned int theTrackID = theTrack->GetTrackID();
158 
159  LogDebug("PLTSensitiveDetector") << " OLD (d,t) = (" << lastId << "," << lastTrack
160  << "), new = (" << theDetUnitId << "," << theTrackID << ") return "
161  << ((theTrackID == lastTrack) && (lastId == theDetUnitId));
162  if ((mySimHit != 0) && (theTrackID == lastTrack) && (lastId == theDetUnitId) && closeHit(aStep))
163  return false;
164  return true;
165 }
#define LogDebug(id)
UpdatablePSimHit * mySimHit
virtual bool closeHit(G4Step *)
virtual uint32_t setDetUnitId(G4Step *)
bool PLTSensitiveDetector::ProcessHits ( G4Step *  aStep,
G4TouchableHistory *   
)
virtual

Implements SensitiveDetector.

Definition at line 69 of file PLTSensitiveDetector.cc.

References createHit(), energyCut, energyHistoryCut, getOrCreateTrackInformation(), info, lastTrack, LogDebug, newHit(), TrackInformation::putInHistory(), sendHit(), TrackInformation::storeTrack(), and updateHit().

69  {
70 
71  LogDebug("PLTSensitiveDetector") << " Entering a new Step "
72  << aStep->GetTotalEnergyDeposit() << " "
73  << aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName();
74 
75  G4Track * gTrack = aStep->GetTrack();
76  if ((unsigned int)(gTrack->GetTrackID()) != lastTrack) {
77 
78  if (gTrack->GetKineticEnergy() > energyCut){
80  info->storeTrack(true);
81  }
82  if (gTrack->GetKineticEnergy() > energyHistoryCut){
84  info->putInHistory();
85  }
86  }
87 
88  if (aStep->GetTotalEnergyDeposit()>0.) {
89  if (newHit(aStep) == true) {
90  sendHit();
91  createHit(aStep);
92  } else {
93  updateHit(aStep);
94  }
95  return true;
96  }
97  return false;
98 }
#define LogDebug(id)
virtual void createHit(G4Step *)
bool storeTrack() const
TrackInformation * getOrCreateTrackInformation(const G4Track *)
virtual bool newHit(G4Step *)
virtual void updateHit(G4Step *)
void PLTSensitiveDetector::sendHit ( )
privatevirtual

Definition at line 121 of file PLTSensitiveDetector.cc.

References PSimHit::detUnitId(), PSimHit::energyLoss(), PSimHit::entryPoint(), PSimHit::exitPoint(), lastId, lastTrack, LogDebug, mySimHit, pname, TrackingSlaveSD::processHits(), slave, and PSimHit::trackId().

Referenced by EndOfEvent(), and ProcessHits().

121  {
122  if (mySimHit == 0) return;
123  LogDebug("PLTSensitiveDetector") << " Storing PSimHit: " << pname << " " << mySimHit->detUnitId()
124  << " " << mySimHit->trackId() << " " << mySimHit->energyLoss()
125  << " " << mySimHit->entryPoint() << " " << mySimHit->exitPoint();
126 
128 
129  // clean up
130  delete mySimHit;
131  mySimHit = 0;
132  lastTrack = 0;
133  lastId = 0;
134 }
#define LogDebug(id)
UpdatablePSimHit * mySimHit
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:38
TrackingSlaveSD * slave
float energyLoss() const
The energy deposit in the PSimHit, in ???.
Definition: PSimHit.h:75
unsigned int trackId() const
Definition: PSimHit.h:102
virtual bool processHits(const PSimHit &)
Local3DPoint entryPoint() const
Entry point in the local Det frame.
Definition: PSimHit.h:35
unsigned int detUnitId() const
Definition: PSimHit.h:93
uint32_t PLTSensitiveDetector::setDetUnitId ( G4Step *  )
virtual

Implements SensitiveDetector.

Definition at line 100 of file PLTSensitiveDetector.cc.

References LogDebug.

Referenced by createHit(), and newHit().

100  {
101 
102  unsigned int detId = 0;
103 
104  LogDebug("PLTSensitiveDetector")<< " DetID = "<<detId;
105 
106  return detId;
107 }
#define LogDebug(id)
void PLTSensitiveDetector::update ( const BeginOfEvent )
privatevirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfEvent * >.

Definition at line 227 of file PLTSensitiveDetector.cc.

References clearHits(), eventno, and mySimHit.

Referenced by progressbar.ProgressBar::__next__(), relval_steps.Matrix::__setitem__(), relval_steps.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(), relval_steps.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().

227  {
228 
229  clearHits();
230  eventno = (*i)()->GetEventID();
231  mySimHit = 0;
232 }
UpdatablePSimHit * mySimHit
void PLTSensitiveDetector::update ( const BeginOfTrack )
privatevirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfTrack * >.

Definition at line 234 of file PLTSensitiveDetector.cc.

References pname.

Referenced by progressbar.ProgressBar::__next__(), relval_steps.Matrix::__setitem__(), relval_steps.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(), relval_steps.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().

234  {
235 
236  const G4Track* gTrack = (*bot)();
237  pname = gTrack->GetDynamicParticle()->GetDefinition()->GetParticleName();
238 }
void PLTSensitiveDetector::update ( const BeginOfJob )
privatevirtual

This routine will be called when the appropriate signal arrives.

Implements Observer< const BeginOfJob * >.

Definition at line 225 of file PLTSensitiveDetector.cc.

Referenced by progressbar.ProgressBar::__next__(), relval_steps.Matrix::__setitem__(), relval_steps.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(), relval_steps.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().

225 { }
void PLTSensitiveDetector::updateHit ( G4Step *  aStep)
privatevirtual

Definition at line 136 of file PLTSensitiveDetector.cc.

References UpdatablePSimHit::addEnergyLoss(), PSimHit::detUnitId(), PSimHit::energyLoss(), PSimHit::entryPoint(), PSimHit::exitPoint(), SensitiveDetector::FinalStepPosition(), globalExitPoint, SensitiveDetector::LocalCoordinates, LogDebug, mySimHit, pname, UpdatablePSimHit::setExitPoint(), PSimHit::trackId(), and SensitiveDetector::WorldCoordinates.

Referenced by ProcessHits().

136  {
137 
139  float theEnergyLoss = aStep->GetTotalEnergyDeposit()/GeV;
140  mySimHit->setExitPoint(theExitPoint);
141  LogDebug("PLTSensitiveDetector")<< " Before : " << mySimHit->energyLoss();
142  mySimHit->addEnergyLoss(theEnergyLoss);
144 
145  LogDebug("PLTSensitiveDetector") << " Updating: new exitpoint " << pname << " "
146  << theExitPoint << " new energy loss " << theEnergyLoss
147  << "\n Updated PSimHit: " << mySimHit->detUnitId()
148  << " " << mySimHit->trackId()
149  << " " << mySimHit->energyLoss() << " "
150  << mySimHit->entryPoint() << " " << mySimHit->exitPoint();
151 }
#define LogDebug(id)
UpdatablePSimHit * mySimHit
void setExitPoint(const Local3DPoint &exit)
Local3DPoint exitPoint() const
Exit point in the local Det frame.
Definition: PSimHit.h:38
void addEnergyLoss(float eloss)
float energyLoss() const
The energy deposit in the PSimHit, in ???.
Definition: PSimHit.h:75
unsigned int trackId() const
Definition: PSimHit.h:102
Local3DPoint FinalStepPosition(G4Step *s, coordinates)
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

float PLTSensitiveDetector::energyCut
private

Definition at line 68 of file PLTSensitiveDetector.h.

Referenced by PLTSensitiveDetector(), and ProcessHits().

float PLTSensitiveDetector::energyHistoryCut
private

Definition at line 69 of file PLTSensitiveDetector.h.

Referenced by PLTSensitiveDetector(), and ProcessHits().

int PLTSensitiveDetector::eventno
private

Definition at line 76 of file PLTSensitiveDetector.h.

Referenced by update().

Local3DPoint PLTSensitiveDetector::globalEntryPoint
private

Definition at line 71 of file PLTSensitiveDetector.h.

Referenced by createHit().

Local3DPoint PLTSensitiveDetector::globalExitPoint
private

Definition at line 72 of file PLTSensitiveDetector.h.

Referenced by createHit(), and updateHit().

uint32_t PLTSensitiveDetector::lastId
private

Definition at line 74 of file PLTSensitiveDetector.h.

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

unsigned int PLTSensitiveDetector::lastTrack
private

Definition at line 75 of file PLTSensitiveDetector.h.

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

G4TrackToParticleID* PLTSensitiveDetector::myG4TrackToParticleID
private

Definition at line 65 of file PLTSensitiveDetector.h.

Referenced by createHit(), PLTSensitiveDetector(), and ~PLTSensitiveDetector().

std::string PLTSensitiveDetector::myName
private

Definition at line 66 of file PLTSensitiveDetector.h.

Referenced by EndOfEvent().

UpdatablePSimHit* PLTSensitiveDetector::mySimHit
private

Definition at line 67 of file PLTSensitiveDetector.h.

Referenced by closeHit(), createHit(), EndOfEvent(), newHit(), sendHit(), update(), and updateHit().

G4VPhysicalVolume* PLTSensitiveDetector::oldVolume
private

Definition at line 73 of file PLTSensitiveDetector.h.

Referenced by createHit().

std::string PLTSensitiveDetector::pname
private

Definition at line 77 of file PLTSensitiveDetector.h.

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

TrackingSlaveSD* PLTSensitiveDetector::slave
private
G4ProcessTypeEnumerator* PLTSensitiveDetector::theG4ProcessTypeEnumerator
private

Definition at line 64 of file PLTSensitiveDetector.h.

Referenced by createHit(), PLTSensitiveDetector(), and ~PLTSensitiveDetector().