#include <PLTSensitiveDetector.h>
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 () |
Private Member Functions | |
virtual void | clearHits () |
virtual bool | closeHit (G4Step *) |
virtual void | createHit (G4Step *) |
TrackInformation * | getOrCreateTrackInformation (const G4Track *) |
virtual bool | newHit (G4Step *) |
virtual void | sendHit () |
void | update (const BeginOfJob *) |
This routine will be called when the appropriate signal arrives. | |
void | update (const BeginOfTrack *) |
This routine will be called when the appropriate signal arrives. | |
void | update (const BeginOfEvent *) |
This routine will be called when the appropriate signal arrives. | |
virtual void | updateHit (G4Step *) |
Private Attributes | |
float | energyCut |
float | energyHistoryCut |
int | eventno |
Local3DPoint | globalEntryPoint |
Local3DPoint | globalExitPoint |
uint32_t | lastId |
unsigned int | lastTrack |
G4TrackToParticleID * | myG4TrackToParticleID |
std::string | myName |
UpdatablePSimHit * | mySimHit |
G4VPhysicalVolume * | oldVolume |
std::string | pname |
TrackingSlaveSD * | slave |
G4ProcessTypeEnumerator * | theG4ProcessTypeEnumerator |
Definition at line 30 of file PLTSensitiveDetector.h.
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.
: SensitiveTkDetector(name, cpv, clg, p), myName(name), mySimHit(0), oldVolume(0), lastId(0), lastTrack(0), eventno(0) { edm::ParameterSet m_TrackerSD = p.getParameter<edm::ParameterSet>("PltSD"); energyCut = m_TrackerSD.getParameter<double>("EnergyThresholdForPersistencyInGeV")*GeV; //default must be 0.5 energyHistoryCut = m_TrackerSD.getParameter<double>("EnergyThresholdForHistoryInGeV")*GeV;//default must be 0.05 edm::LogInfo("PLTSensitiveDetector") <<"Criteria for Saving Tracker SimTracks: \n " <<" History: "<<energyHistoryCut<< " MeV ; Persistency: "<< energyCut<<" MeV\n" <<" Constructing a PLTSensitiveDetector with "; slave = new TrackingSlaveSD(name); // Now attach the right detectors (LogicalVolumes) to me std::vector<std::string> lvNames = clg.logicalNames(name); this->Register(); for (std::vector<std::string>::iterator it = lvNames.begin(); it != lvNames.end(); it++) { edm::LogInfo("PLTSensitiveDetector")<< name << " attaching LV " << *it; this->AssignSD(*it); } theG4ProcessTypeEnumerator = new G4ProcessTypeEnumerator; myG4TrackToParticleID = new G4TrackToParticleID; }
PLTSensitiveDetector::~PLTSensitiveDetector | ( | ) | [virtual] |
Definition at line 63 of file PLTSensitiveDetector.cc.
References myG4TrackToParticleID, slave, and theG4ProcessTypeEnumerator.
{ delete slave; delete theG4ProcessTypeEnumerator; delete myG4TrackToParticleID; }
void PLTSensitiveDetector::clearHits | ( | ) | [private, virtual] |
Implements SensitiveDetector.
Definition at line 240 of file PLTSensitiveDetector.cc.
References TrackingSlaveSD::Initialize(), and slave.
Referenced by update().
{ slave->Initialize(); }
bool PLTSensitiveDetector::closeHit | ( | G4Step * | aStep | ) | [private, virtual] |
Definition at line 167 of file PLTSensitiveDetector.cc.
References PSimHit::exitPoint(), SensitiveDetector::InitialStepPosition(), SensitiveDetector::LocalCoordinates, LogDebug, mag(), and mySimHit.
Referenced by newHit().
{ if (mySimHit == 0) return false; const float tolerance = 0.05 * mm; // 50 micron are allowed between the exit // point of the current hit and the entry point of the new hit Local3DPoint theEntryPoint = SensitiveDetector::InitialStepPosition(aStep,LocalCoordinates); LogDebug("PLTSensitiveDetector")<< " closeHit: distance = " << (mySimHit->exitPoint()-theEntryPoint).mag(); if ((mySimHit->exitPoint()-theEntryPoint).mag()<tolerance) return true; return false; }
void PLTSensitiveDetector::createHit | ( | G4Step * | aStep | ) | [private, virtual] |
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(), v, and SensitiveDetector::WorldCoordinates.
Referenced by ProcessHits().
{ if (mySimHit != 0) { delete mySimHit; mySimHit=0; } G4Track * theTrack = aStep->GetTrack(); G4VPhysicalVolume * v = aStep->GetPreStepPoint()->GetPhysicalVolume(); Local3DPoint theEntryPoint = SensitiveDetector::InitialStepPosition(aStep,LocalCoordinates); Local3DPoint theExitPoint = SensitiveDetector::FinalStepPosition(aStep,LocalCoordinates); float thePabs = aStep->GetPreStepPoint()->GetMomentum().mag()/GeV; float theTof = aStep->GetPreStepPoint()->GetGlobalTime()/nanosecond; float theEnergyLoss = aStep->GetTotalEnergyDeposit()/GeV; int theParticleType = myG4TrackToParticleID->particleID(theTrack); uint32_t theDetUnitId = setDetUnitId(aStep); globalEntryPoint = SensitiveDetector::InitialStepPosition(aStep,WorldCoordinates); globalExitPoint = SensitiveDetector::FinalStepPosition(aStep,WorldCoordinates); pname = theTrack->GetDynamicParticle()->GetDefinition()->GetParticleName(); unsigned int theTrackID = theTrack->GetTrackID(); G4ThreeVector gmd = aStep->GetPreStepPoint()->GetMomentumDirection(); // convert it to local frame G4ThreeVector lmd = ((G4TouchableHistory *)(aStep->GetPreStepPoint()->GetTouchable()))->GetHistory()->GetTopTransform().TransformAxis(gmd); Local3DPoint lnmd = ConvertToLocal3DPoint(lmd); float theThetaAtEntry = lnmd.theta(); float thePhiAtEntry = lnmd.phi(); mySimHit = new UpdatablePSimHit(theEntryPoint,theExitPoint,thePabs,theTof, theEnergyLoss,theParticleType,theDetUnitId, theTrackID,theThetaAtEntry,thePhiAtEntry, theG4ProcessTypeEnumerator->processId(theTrack->GetCreatorProcess())); LogDebug("PLTSensitiveDetector") << " Created PSimHit: " << pname << " " << mySimHit->detUnitId() << " " << mySimHit->trackId() << " " << mySimHit->energyLoss() << " " << mySimHit->entryPoint() << " " << mySimHit->exitPoint(); lastId = theDetUnitId; lastTrack = theTrackID; oldVolume = v; }
void PLTSensitiveDetector::EndOfEvent | ( | G4HCofThisEvent * | ) | [virtual] |
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.
TrackInformation * PLTSensitiveDetector::getOrCreateTrackInformation | ( | const G4Track * | gTrack | ) | [private] |
Definition at line 244 of file PLTSensitiveDetector.cc.
References info, and groupFilesInBlocks::temp.
Referenced by ProcessHits().
{ G4VUserTrackInformation* temp = gTrack->GetUserInformation(); if (temp == 0){ edm::LogError("PLTSensitiveDetector") <<" ERROR: no G4VUserTrackInformation available"; abort(); }else{ TrackInformation* info = dynamic_cast<TrackInformation*>(temp); if (info == 0){ edm::LogError("PLTSensitiveDetector") <<" ERROR: TkSimTrackSelection: the UserInformation does not appear to be a TrackInformation"; abort(); } return info; } }
bool PLTSensitiveDetector::newHit | ( | G4Step * | aStep | ) | [private, virtual] |
Definition at line 153 of file PLTSensitiveDetector.cc.
References closeHit(), lastId, lastTrack, LogDebug, mySimHit, and setDetUnitId().
Referenced by ProcessHits().
{ G4Track * theTrack = aStep->GetTrack(); uint32_t theDetUnitId = setDetUnitId(aStep); unsigned int theTrackID = theTrack->GetTrackID(); LogDebug("PLTSensitiveDetector") << " OLD (d,t) = (" << lastId << "," << lastTrack << "), new = (" << theDetUnitId << "," << theTrackID << ") return " << ((theTrackID == lastTrack) && (lastId == theDetUnitId)); if ((mySimHit != 0) && (theTrackID == lastTrack) && (lastId == theDetUnitId) && closeHit(aStep)) return false; return true; }
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().
{ LogDebug("PLTSensitiveDetector") << " Entering a new Step " << aStep->GetTotalEnergyDeposit() << " " << aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName(); G4Track * gTrack = aStep->GetTrack(); if ((unsigned int)(gTrack->GetTrackID()) != lastTrack) { if (gTrack->GetKineticEnergy() > energyCut){ TrackInformation* info = getOrCreateTrackInformation(gTrack); info->storeTrack(true); } if (gTrack->GetKineticEnergy() > energyHistoryCut){ TrackInformation* info = getOrCreateTrackInformation(gTrack); info->putInHistory(); } } if (aStep->GetTotalEnergyDeposit()>0.) { if (newHit(aStep) == true) { sendHit(); createHit(aStep); } else { updateHit(aStep); } return true; } return false; }
void PLTSensitiveDetector::sendHit | ( | ) | [private, virtual] |
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().
{ if (mySimHit == 0) return; LogDebug("PLTSensitiveDetector") << " Storing PSimHit: " << pname << " " << mySimHit->detUnitId() << " " << mySimHit->trackId() << " " << mySimHit->energyLoss() << " " << mySimHit->entryPoint() << " " << mySimHit->exitPoint(); slave->processHits(*mySimHit); // clean up delete mySimHit; mySimHit = 0; lastTrack = 0; lastId = 0; }
uint32_t PLTSensitiveDetector::setDetUnitId | ( | G4Step * | ) | [virtual] |
Implements SensitiveDetector.
Definition at line 100 of file PLTSensitiveDetector.cc.
References LogDebug.
Referenced by createHit(), and newHit().
{ unsigned int detId = 0; LogDebug("PLTSensitiveDetector")<< " DetID = "<<detId; return detId; }
void PLTSensitiveDetector::update | ( | const BeginOfTrack * | ) | [private, virtual] |
This routine will be called when the appropriate signal arrives.
Implements Observer< const BeginOfTrack * >.
Definition at line 234 of file PLTSensitiveDetector.cc.
References pname.
{ const G4Track* gTrack = (*bot)(); pname = gTrack->GetDynamicParticle()->GetDefinition()->GetParticleName(); }
void PLTSensitiveDetector::update | ( | const BeginOfJob * | ) | [private, virtual] |
This routine will be called when the appropriate signal arrives.
Implements Observer< const BeginOfJob * >.
Definition at line 225 of file PLTSensitiveDetector.cc.
{ }
void PLTSensitiveDetector::update | ( | const BeginOfEvent * | ) | [private, virtual] |
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.
void PLTSensitiveDetector::updateHit | ( | G4Step * | aStep | ) | [private, virtual] |
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().
{ Local3DPoint theExitPoint = SensitiveDetector::FinalStepPosition(aStep,LocalCoordinates); float theEnergyLoss = aStep->GetTotalEnergyDeposit()/GeV; mySimHit->setExitPoint(theExitPoint); LogDebug("PLTSensitiveDetector")<< " Before : " << mySimHit->energyLoss(); mySimHit->addEnergyLoss(theEnergyLoss); globalExitPoint = SensitiveDetector::FinalStepPosition(aStep,WorldCoordinates); LogDebug("PLTSensitiveDetector") << " Updating: new exitpoint " << pname << " " << theExitPoint << " new energy loss " << theEnergyLoss << "\n Updated PSimHit: " << mySimHit->detUnitId() << " " << mySimHit->trackId() << " " << mySimHit->energyLoss() << " " << mySimHit->entryPoint() << " " << mySimHit->exitPoint(); }
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().
Definition at line 71 of file PLTSensitiveDetector.h.
Referenced by createHit().
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().
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] |
Definition at line 63 of file PLTSensitiveDetector.h.
Referenced by clearHits(), fillHits(), PLTSensitiveDetector(), sendHit(), and ~PLTSensitiveDetector().
Definition at line 64 of file PLTSensitiveDetector.h.
Referenced by createHit(), PLTSensitiveDetector(), and ~PLTSensitiveDetector().