#include <TkAccumulatingSensitiveDetector.h>
Public Member Functions | |
virtual void | EndOfEvent (G4HCofThisEvent *) |
void | fillHits (edm::PSimHitContainer &, std::string use) |
std::vector< std::string > | getNames () |
virtual bool | ProcessHits (G4Step *, G4TouchableHistory *) |
virtual uint32_t | setDetUnitId (G4Step *) |
TkAccumulatingSensitiveDetector (std::string, const DDCompactView &, SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *) | |
std::string | type () |
virtual | ~TkAccumulatingSensitiveDetector () |
Private Member Functions | |
void | checkExitPoint (Local3DPoint) |
virtual void | clearHits () |
virtual bool | closeHit (G4Step *) |
virtual void | createHit (G4Step *) |
TrackInformation * | getOrCreateTrackInformation (const G4Track *) |
virtual bool | newHit (G4Step *) |
virtual void | sendHit () |
int | tofBin (float) |
Local3DPoint | toOrcaRef (Local3DPoint, G4VPhysicalVolume *) |
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 | |
bool | allowZeroEnergyLoss |
float | energyCut |
float | energyHistoryCut |
int | eventno |
Local3DPoint | globalEntryPoint |
Local3DPoint | globalExitPoint |
uint32_t | lastId |
unsigned int | lastTrack |
G4TrackToParticleID * | myG4TrackToParticleID |
std::string | myName |
FrameRotation * | myRotation |
UpdatablePSimHit * | mySimHit |
bool | neverAccumulate |
TrackerG4SimHitNumberingScheme * | numberingScheme_ |
G4VPhysicalVolume * | oldVolume |
std::string | pname |
bool | printHits |
float | px |
float | py |
float | pz |
float | rTracker |
TrackingSlaveSD * | slaveHighTof |
TrackingSlaveSD * | slaveLowTof |
G4ProcessTypeEnumerator * | theG4ProcessTypeEnumerator |
const SimTrackManager * | theManager |
double | theSigma |
float | zTracker |
Definition at line 30 of file TkAccumulatingSensitiveDetector.h.
TkAccumulatingSensitiveDetector::TkAccumulatingSensitiveDetector | ( | std::string | , |
const DDCompactView & | , | ||
SensitiveDetectorCatalog & | , | ||
edm::ParameterSet const & | , | ||
const SimTrackManager * | |||
) |
Definition at line 61 of file TkAccumulatingSensitiveDetector.cc.
References allowZeroEnergyLoss, SensitiveDetector::AssignSD(), energyCut, energyHistoryCut, edm::ParameterSet::getParameter(), SensitiveDetectorCatalog::logicalNames(), myG4TrackToParticleID, myName, myRotation, neverAccumulate, printHits, SensitiveDetector::Register(), slaveHighTof, slaveLowTof, theG4ProcessTypeEnumerator, and theSigma.
: SensitiveTkDetector(name, cpv, clg, p), myName(name), myRotation(0), mySimHit(0),theManager(manager), oldVolume(0), lastId(0), lastTrack(0), eventno(0) ,rTracker(1200.*mm),zTracker(3000.*mm), numberingScheme_(0) { edm::ParameterSet m_TrackerSD = p.getParameter<edm::ParameterSet>("TrackerSD"); allowZeroEnergyLoss = m_TrackerSD.getParameter<bool>("ZeroEnergyLoss"); neverAccumulate = m_TrackerSD.getParameter<bool>("NeverAccumulate"); printHits = m_TrackerSD.getParameter<bool>("PrintHits"); theSigma = m_TrackerSD.getParameter<double>("ElectronicSigmaInNanoSeconds"); 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("TrackerSimInfo") <<"Criteria for Saving Tracker SimTracks: "; edm::LogInfo("TrackerSimInfo")<<" History: "<<energyHistoryCut<< " MeV ; Persistency: "<< energyCut<<" MeV "; edm::LogInfo("TrackerSimInfo")<<" Constructing a TkAccumulatingSensitiveDetector with "; #ifndef FAKEFRAMEROTATION // No Rotation given in input, automagically choose one based upon the name if (name.find("TrackerHits") != string::npos) { edm::LogInfo("TrackerSimInfo")<<" TkAccumulatingSensitiveDetector: using TrackerFrameRotation for "<<myName; myRotation = new TrackerFrameRotation; } // Just in case (test beam etc) if (myRotation == 0) { edm::LogInfo("TrackerSimInfo")<<" TkAccumulatingSensitiveDetector: using StandardFrameRotation for "<<myName; myRotation = new StandardFrameRotation; } #else myRotation = new FakeFrameRotation; edm::LogWarning("TrackerSimInfo")<< " WARNING - Using FakeFrameRotation in TkAccumulatingSensitiveDetector;"; #endif slaveLowTof = new TrackingSlaveSD(name+"LowTof"); slaveHighTof = new TrackingSlaveSD(name+"HighTof"); // Now attach the right detectors (LogicalVolumes) to me vector<string> lvNames = clg.logicalNames(name); this->Register(); for (vector<string>::iterator it = lvNames.begin(); it != lvNames.end(); it++) { edm::LogInfo("TrackerSimInfo")<< name << " attaching LV " << *it; this->AssignSD(*it); } theG4ProcessTypeEnumerator = new G4ProcessTypeEnumerator; myG4TrackToParticleID = new G4TrackToParticleID; }
TkAccumulatingSensitiveDetector::~TkAccumulatingSensitiveDetector | ( | ) | [virtual] |
Definition at line 117 of file TkAccumulatingSensitiveDetector.cc.
References myG4TrackToParticleID, slaveHighTof, slaveLowTof, and theG4ProcessTypeEnumerator.
{ delete slaveLowTof; delete slaveHighTof; delete theG4ProcessTypeEnumerator; delete myG4TrackToParticleID; }
void TkAccumulatingSensitiveDetector::checkExitPoint | ( | Local3DPoint | p | ) | [private] |
Definition at line 424 of file TkAccumulatingSensitiveDetector.cc.
References abs, PV3DBase< T, PVType, FrameType >::z(), and z.
Referenced by createHit(), and updateHit().
{ double z = p.z(); if (std::abs(z)<0.3*mm) return; bool sendExc = false; //static SimpleConfigurable<bool> sendExc(false,"TrackerSim:ThrowOnBadHits"); edm::LogWarning("TrackerSimInfo")<< " ************ Hit outside the detector ; Local z " << z << "; skipping event = " << sendExc; if (sendExc == true) { G4EventManager::GetEventManager()->AbortCurrentEvent(); G4EventManager::GetEventManager()->GetNonconstCurrentEvent()->SetEventAborted(); } }
void TkAccumulatingSensitiveDetector::clearHits | ( | ) | [private, virtual] |
Implements SensitiveDetector.
Definition at line 418 of file TkAccumulatingSensitiveDetector.cc.
References TrackingSlaveSD::Initialize(), slaveHighTof, and slaveLowTof.
Referenced by update().
{ slaveLowTof->Initialize(); slaveHighTof->Initialize(); }
bool TkAccumulatingSensitiveDetector::closeHit | ( | G4Step * | aStep | ) | [private, virtual] |
Definition at line 377 of file TkAccumulatingSensitiveDetector.cc.
References PSimHit::exitPoint(), SensitiveDetector::InitialStepPosition(), SensitiveDetector::LocalCoordinates, LogDebug, mag(), mySimHit, toOrcaRef(), and v.
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 G4VPhysicalVolume * v = aStep->GetPreStepPoint()->GetPhysicalVolume(); Local3DPoint theEntryPoint = toOrcaRef(SensitiveDetector::InitialStepPosition(aStep,LocalCoordinates),v); LogDebug("TrackerSimDebug")<< " closeHit: distance = " << (mySimHit->exitPoint()-theEntryPoint).mag(); if ((mySimHit->exitPoint()-theEntryPoint).mag()<tolerance) return true; return false; }
void TkAccumulatingSensitiveDetector::createHit | ( | G4Step * | aStep | ) | [private, virtual] |
Definition at line 255 of file TkAccumulatingSensitiveDetector.cc.
References checkExitPoint(), SensitiveDetector::ConvertToLocal3DPoint(), PSimHit::detUnitId(), PSimHit::energyLoss(), PSimHit::entryPoint(), PSimHit::exitPoint(), SensitiveDetector::FinalStepPosition(), globalEntryPoint, globalExitPoint, info, SensitiveDetector::InitialStepPosition(), lastId, lastTrack, SensitiveDetector::LocalCoordinates, LogDebug, myG4TrackToParticleID, mySimHit, oldVolume, G4TrackToParticleID::particleID(), PV3DBase< T, PVType, FrameType >::phi(), pname, G4ProcessTypeEnumerator::processId(), px, py, pz, setDetUnitId(), TrackInformation::storeTrack(), groupFilesInBlocks::temp, theG4ProcessTypeEnumerator, PV3DBase< T, PVType, FrameType >::theta(), toOrcaRef(), 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 = toOrcaRef(SensitiveDetector::InitialStepPosition(aStep,LocalCoordinates),v); Local3DPoint theExitPoint = toOrcaRef(SensitiveDetector::FinalStepPosition(aStep,LocalCoordinates),v); // // This allows to send he skipEvent if it is outside! // checkExitPoint(theExitPoint); 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(); if (theDetUnitId == 0) { edm::LogError("TrackerSimInfo") << " Error: theDetUnitId is not valid."; abort(); } unsigned int theTrackID = theTrack->GetTrackID(); // To whom assign the Hit? // First iteration: if the track is to be stored, use the current number; // otherwise, get to the mother unsigned int theTrackIDInsideTheSimHit=theTrackID; G4VUserTrackInformation * info = theTrack->GetUserInformation(); if (info == 0) edm::LogError("TrackerSimInfo")<< " Error: no UserInformation available "; else { TrackInformation * temp = dynamic_cast<TrackInformation* >(info); if (temp ==0) edm::LogError("TrackerSimInfo")<< " Error:G4VUserTrackInformation is not a TrackInformation."; if (temp->storeTrack() == false) { // Go to the mother! LogDebug("TrackerSimDebug")<< " TkAccumulatingSensitiveDetector:createHit(): setting the TrackID from " << theTrackIDInsideTheSimHit; theTrackIDInsideTheSimHit = theTrack->GetParentID(); LogDebug("TrackerSimDebug")<< " to the mother one " << theTrackIDInsideTheSimHit << " " << theEnergyLoss; } else { LogDebug("TrackerSimDebug")<< " TkAccumulatingSensitiveDetector:createHit(): leaving the current TrackID " << theTrackIDInsideTheSimHit; } } px = aStep->GetPreStepPoint()->GetMomentum().x()/GeV; py = aStep->GetPreStepPoint()->GetMomentum().y()/GeV; pz = aStep->GetPreStepPoint()->GetMomentum().z()/GeV; G4ThreeVector gmd = aStep->GetPreStepPoint()->GetMomentumDirection(); // convert it to local frame G4ThreeVector lmd = ((G4TouchableHistory *)(aStep->GetPreStepPoint()->GetTouchable()))->GetHistory() ->GetTopTransform().TransformAxis(gmd); Local3DPoint lnmd = toOrcaRef(ConvertToLocal3DPoint(lmd),v); float theThetaAtEntry = lnmd.theta(); float thePhiAtEntry = lnmd.phi(); mySimHit = new UpdatablePSimHit(theEntryPoint,theExitPoint,thePabs,theTof, theEnergyLoss,theParticleType,theDetUnitId, theTrackIDInsideTheSimHit,theThetaAtEntry,thePhiAtEntry, theG4ProcessTypeEnumerator->processId(theTrack->GetCreatorProcess())); LogDebug("TrackerSimDebug")<< " Created PSimHit: " << pname << " " << mySimHit->detUnitId() << " " << mySimHit->trackId() << " " << mySimHit->energyLoss() << " " << mySimHit->entryPoint() << " " << mySimHit->exitPoint(); lastId = theDetUnitId; lastTrack = theTrackID; oldVolume = v; }
void TkAccumulatingSensitiveDetector::EndOfEvent | ( | G4HCofThisEvent * | ) | [virtual] |
void TkAccumulatingSensitiveDetector::fillHits | ( | edm::PSimHitContainer & | c, |
std::string | use | ||
) | [virtual] |
Implements SensitiveTkDetector.
Definition at line 454 of file TkAccumulatingSensitiveDetector.cc.
References TrackingSlaveSD::hits(), n, TrackingSlaveSD::name(), slaveHighTof, and slaveLowTof.
{ // // do it once for low, once for High // if (slaveLowTof->name() == n) c=slaveLowTof->hits(); if (slaveHighTof->name() == n) c=slaveHighTof->hits(); }
std::vector< string > TkAccumulatingSensitiveDetector::getNames | ( | ) | [virtual] |
Reimplemented from SensitiveDetector.
Definition at line 464 of file TkAccumulatingSensitiveDetector.cc.
References TrackingSlaveSD::name(), slaveHighTof, slaveLowTof, and groupFilesInBlocks::temp.
{ std::vector<string> temp; temp.push_back(slaveLowTof->name()); temp.push_back(slaveHighTof->name()); return temp; }
TrackInformation * TkAccumulatingSensitiveDetector::getOrCreateTrackInformation | ( | const G4Track * | gTrack | ) | [private] |
Definition at line 439 of file TkAccumulatingSensitiveDetector.cc.
References info, and groupFilesInBlocks::temp.
Referenced by update().
{ G4VUserTrackInformation* temp = gTrack->GetUserInformation(); if (temp == 0){ edm::LogError("TrackerSimInfo") <<" ERROR: no G4VUserTrackInformation available"; abort(); }else{ TrackInformation* info = dynamic_cast<TrackInformation*>(temp); if (info ==0){ edm::LogError("TrackerSimInfo")<<" ERROR: TkSimTrackSelection: the UserInformation does not appear to be a TrackInformation"; abort(); } return info; } }
bool TkAccumulatingSensitiveDetector::newHit | ( | G4Step * | aStep | ) | [private, virtual] |
Definition at line 362 of file TkAccumulatingSensitiveDetector.cc.
References closeHit(), lastId, lastTrack, LogDebug, mySimHit, neverAccumulate, and setDetUnitId().
Referenced by ProcessHits().
{ if (neverAccumulate == true) return true; G4Track * theTrack = aStep->GetTrack(); uint32_t theDetUnitId = setDetUnitId(aStep); unsigned int theTrackID = theTrack->GetTrackID(); LogDebug("TrackerSimDebug")<< " 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 TkAccumulatingSensitiveDetector::ProcessHits | ( | G4Step * | aStep, |
G4TouchableHistory * | ROhist | ||
) | [virtual] |
Implements SensitiveDetector.
Definition at line 126 of file TkAccumulatingSensitiveDetector.cc.
References allowZeroEnergyLoss, createHit(), LogDebug, newHit(), sendHit(), and updateHit().
Referenced by LaserAlignmentSimulation::update().
{ LogDebug("TrackerSimDebug")<< " Entering a new Step " << aStep->GetTotalEnergyDeposit() << " " << aStep->GetPreStepPoint()->GetPhysicalVolume()->GetLogicalVolume()->GetName(); //TimeMe t1(theTimer, false); if (aStep->GetTotalEnergyDeposit()>0. || allowZeroEnergyLoss == true) { if (newHit(aStep) == true) { sendHit(); createHit(aStep); } else { updateHit(aStep); } return true; } return false; }
void TkAccumulatingSensitiveDetector::sendHit | ( | ) | [private, virtual] |
Definition at line 224 of file TkAccumulatingSensitiveDetector.cc.
References PSimHit::detUnitId(), PSimHit::energyLoss(), PSimHit::entryPoint(), eventno, PSimHit::exitPoint(), TkSimHitPrinter::getPropagationSign(), globalEntryPoint, globalExitPoint, lastId, lastTrack, LogDebug, myName, mySimHit, oldVolume, PSimHit::pabs(), pname, TkSimHitPrinter::printGlobal(), TkSimHitPrinter::printGlobalMomentum(), TkSimHitPrinter::printHitData(), printHits, TkSimHitPrinter::printLocal(), TkSimHitPrinter::printMomentumOfTrack(), TrackingSlaveSD::processHits(), px, py, pz, slaveHighTof, slaveLowTof, TkSimHitPrinter::startNewSimHit(), PSimHit::timeOfFlight(), tofBin(), and PSimHit::trackId().
Referenced by EndOfEvent(), and ProcessHits().
{ if (mySimHit == 0) return; LogDebug("TrackerSimDebug")<< " Storing PSimHit: " << pname << " " << mySimHit->detUnitId() << " " << mySimHit->trackId() << " " << mySimHit->energyLoss() << " " << mySimHit->entryPoint() << " " << mySimHit->exitPoint(); if (printHits) { TkSimHitPrinter thePrinter("TkHitPositionOSCAR.dat"); thePrinter.startNewSimHit(myName,oldVolume->GetLogicalVolume()->GetName(), mySimHit->detUnitId(),mySimHit->trackId(),eventno); thePrinter.printLocal(mySimHit->entryPoint(),mySimHit->exitPoint()); thePrinter.printGlobal(globalEntryPoint,globalExitPoint); thePrinter.printHitData(mySimHit->energyLoss(),mySimHit->timeOfFlight()); thePrinter.printMomentumOfTrack(mySimHit->pabs(),pname, thePrinter.getPropagationSign(globalEntryPoint,globalExitPoint)); thePrinter.printGlobalMomentum(px,py,pz); } if (tofBin(mySimHit->timeOfFlight()) == 1) slaveLowTof->processHits(*mySimHit); // implicit conversion (slicing) to PSimHit!!! else slaveHighTof->processHits(*mySimHit); // implicit conversion (slicing) to PSimHit!!! // // clean up delete mySimHit; mySimHit = 0; lastTrack = 0; lastId = 0; }
uint32_t TkAccumulatingSensitiveDetector::setDetUnitId | ( | G4Step * | s | ) | [virtual] |
Implements SensitiveDetector.
Definition at line 149 of file TkAccumulatingSensitiveDetector.cc.
References TrackerG4SimHitNumberingScheme::g4ToNumberingScheme(), LogDebug, and numberingScheme_.
Referenced by createHit(), and newHit().
{ unsigned int detId = numberingScheme_->g4ToNumberingScheme(s->GetPreStepPoint()->GetTouchable()); LogDebug("TrackerSimDebug")<< " DetID = "<<detId; return detId; }
int TkAccumulatingSensitiveDetector::tofBin | ( | float | tof | ) | [private] |
Definition at line 159 of file TkAccumulatingSensitiveDetector.cc.
References theSigma, and dtDQMClient_cfg::threshold.
Referenced by sendHit().
Local3DPoint TkAccumulatingSensitiveDetector::toOrcaRef | ( | Local3DPoint | in, |
G4VPhysicalVolume * | v | ||
) | [private] |
Definition at line 167 of file TkAccumulatingSensitiveDetector.cc.
References myRotation, and FrameRotation::transformPoint().
Referenced by closeHit(), createHit(), and updateHit().
{ if (myRotation !=0) return myRotation->transformPoint(in,v); return (in); }
std::string TkAccumulatingSensitiveDetector::type | ( | ) |
void TkAccumulatingSensitiveDetector::update | ( | const BeginOfEvent * | ) | [private, virtual] |
This routine will be called when the appropriate signal arrives.
Implements Observer< const BeginOfEvent * >.
Definition at line 399 of file TkAccumulatingSensitiveDetector.cc.
References clearHits(), eventno, and mySimHit.
void TkAccumulatingSensitiveDetector::update | ( | const BeginOfJob * | ) | [private, virtual] |
This routine will be called when the appropriate signal arrives.
Implements Observer< const BeginOfJob * >.
Definition at line 406 of file TkAccumulatingSensitiveDetector.cc.
References edm::EventSetup::get(), numberingScheme(), and numberingScheme_.
{ edm::ESHandle<GeometricDet> pDD; const edm::EventSetup* es = (*i)(); es->get<IdealGeometryRecord>().get( pDD ); edm::ESTransientHandle<DDCompactView> pView; es->get<IdealGeometryRecord>().get(pView); numberingScheme_=&(numberingScheme(*pView,*pDD)); }
void TkAccumulatingSensitiveDetector::update | ( | const BeginOfTrack * | ) | [private, virtual] |
This routine will be called when the appropriate signal arrives.
Implements Observer< const BeginOfTrack * >.
Definition at line 174 of file TkAccumulatingSensitiveDetector.cc.
References energyCut, energyHistoryCut, getOrCreateTrackInformation(), info, LogDebug, pos, TrackInformation::putInHistory(), rTracker, TrackInformation::storeTrack(), and zTracker.
{ const G4Track* gTrack = (*bot)(); #ifdef DUMPPROCESSES edm::LogInfo("TrackerSimInfo") <<" -> process creator pointer "<<gTrack->GetCreatorProcess(); if (gTrack->GetCreatorProcess()) edm::LogInfo("TrackerSimInfo")<<" -> PROCESS CREATOR : "<<gTrack->GetCreatorProcess()->GetProcessName(); #endif // //Position // const G4ThreeVector pos = gTrack->GetPosition(); //LogDebug("TrackerSimDebug")<<" ENERGY MeV "<<gTrack->GetKineticEnergy()<<" Energy Cut" << energyCut; //LogDebug("TrackerSimDebug")<<" TOTAL ENERGY "<<gTrack->GetTotalEnergy(); //LogDebug("TrackerSimDebug")<<" WEIGHT "<<gTrack->GetWeight(); // // Check if in Tracker Volume // if (pos.perp() < rTracker &&std::fabs(pos.z()) < zTracker){ // // inside the Tracker // LogDebug("TrackerSimDebug")<<" INSIDE TRACKER"; if (gTrack->GetKineticEnergy() > energyCut){ TrackInformation* info = getOrCreateTrackInformation(gTrack); //LogDebug("TrackerSimDebug")<<" POINTER "<<info; //LogDebug("TrackerSimDebug")<<" track inside the tracker selected for STORE"; //LogDebug("TrackerSimDebug")<<"Track ID (persistent track) = "<<gTrack->GetTrackID(); info->storeTrack(true); } // // Save History? // if (gTrack->GetKineticEnergy() > energyHistoryCut){ TrackInformation* info = getOrCreateTrackInformation(gTrack); info->putInHistory(); //LogDebug("TrackerSimDebug")<<" POINTER "<<info; //LogDebug("TrackerSimDebug")<<" track inside the tracker selected for HISTORY"; //LogDebug("TrackerSimDebug")<<"Track ID (history track) = "<<gTrack->GetTrackID(); } } }
void TkAccumulatingSensitiveDetector::updateHit | ( | G4Step * | aStep | ) | [private, virtual] |
Definition at line 341 of file TkAccumulatingSensitiveDetector.cc.
References UpdatablePSimHit::addEnergyLoss(), checkExitPoint(), PSimHit::detUnitId(), PSimHit::energyLoss(), PSimHit::entryPoint(), PSimHit::exitPoint(), SensitiveDetector::FinalStepPosition(), globalExitPoint, SensitiveDetector::LocalCoordinates, LogDebug, mySimHit, pname, UpdatablePSimHit::setExitPoint(), toOrcaRef(), PSimHit::trackId(), v, and SensitiveDetector::WorldCoordinates.
Referenced by ProcessHits().
{ G4VPhysicalVolume * v = aStep->GetPreStepPoint()->GetPhysicalVolume(); Local3DPoint theExitPoint = toOrcaRef(SensitiveDetector::FinalStepPosition(aStep,LocalCoordinates),v); // // This allows to send he skipEvent if it is outside! // checkExitPoint(theExitPoint); float theEnergyLoss = aStep->GetTotalEnergyDeposit()/GeV; mySimHit->setExitPoint(theExitPoint); LogDebug("TrackerSimDebug")<< " Before : " << mySimHit->energyLoss(); mySimHit->addEnergyLoss(theEnergyLoss); globalExitPoint = SensitiveDetector::FinalStepPosition(aStep,WorldCoordinates); LogDebug("TrackerSimDebug")<< " Updating: new exitpoint " << pname << " " << theExitPoint << " new energy loss " << theEnergyLoss; LogDebug("TrackerSimDebug")<< " Updated PSimHit: " << mySimHit->detUnitId() << " " << mySimHit->trackId() << " " << mySimHit->energyLoss() << " " << mySimHit->entryPoint() << " " << mySimHit->exitPoint(); }
bool TkAccumulatingSensitiveDetector::allowZeroEnergyLoss [private] |
Definition at line 80 of file TkAccumulatingSensitiveDetector.h.
Referenced by ProcessHits(), and TkAccumulatingSensitiveDetector().
float TkAccumulatingSensitiveDetector::energyCut [private] |
Definition at line 85 of file TkAccumulatingSensitiveDetector.h.
Referenced by TkAccumulatingSensitiveDetector(), and update().
float TkAccumulatingSensitiveDetector::energyHistoryCut [private] |
Definition at line 86 of file TkAccumulatingSensitiveDetector.h.
Referenced by TkAccumulatingSensitiveDetector(), and update().
int TkAccumulatingSensitiveDetector::eventno [private] |
Definition at line 77 of file TkAccumulatingSensitiveDetector.h.
Definition at line 69 of file TkAccumulatingSensitiveDetector.h.
Referenced by createHit(), and sendHit().
Definition at line 70 of file TkAccumulatingSensitiveDetector.h.
Referenced by createHit(), sendHit(), and updateHit().
uint32_t TkAccumulatingSensitiveDetector::lastId [private] |
Definition at line 75 of file TkAccumulatingSensitiveDetector.h.
Referenced by createHit(), newHit(), and sendHit().
unsigned int TkAccumulatingSensitiveDetector::lastTrack [private] |
Definition at line 76 of file TkAccumulatingSensitiveDetector.h.
Referenced by createHit(), newHit(), and sendHit().
Definition at line 83 of file TkAccumulatingSensitiveDetector.h.
Referenced by createHit(), TkAccumulatingSensitiveDetector(), and ~TkAccumulatingSensitiveDetector().
std::string TkAccumulatingSensitiveDetector::myName [private] |
Definition at line 63 of file TkAccumulatingSensitiveDetector.h.
Referenced by EndOfEvent(), sendHit(), and TkAccumulatingSensitiveDetector().
Definition at line 66 of file TkAccumulatingSensitiveDetector.h.
Referenced by TkAccumulatingSensitiveDetector(), and toOrcaRef().
Definition at line 67 of file TkAccumulatingSensitiveDetector.h.
Referenced by closeHit(), createHit(), EndOfEvent(), newHit(), sendHit(), update(), and updateHit().
bool TkAccumulatingSensitiveDetector::neverAccumulate [private] |
Definition at line 82 of file TkAccumulatingSensitiveDetector.h.
Referenced by newHit(), and TkAccumulatingSensitiveDetector().
Definition at line 93 of file TkAccumulatingSensitiveDetector.h.
Referenced by setDetUnitId(), and update().
G4VPhysicalVolume* TkAccumulatingSensitiveDetector::oldVolume [private] |
Definition at line 72 of file TkAccumulatingSensitiveDetector.h.
Referenced by createHit(), and sendHit().
std::string TkAccumulatingSensitiveDetector::pname [private] |
Definition at line 68 of file TkAccumulatingSensitiveDetector.h.
Referenced by createHit(), sendHit(), and updateHit().
bool TkAccumulatingSensitiveDetector::printHits [private] |
Definition at line 81 of file TkAccumulatingSensitiveDetector.h.
Referenced by sendHit(), and TkAccumulatingSensitiveDetector().
float TkAccumulatingSensitiveDetector::px [private] |
Definition at line 79 of file TkAccumulatingSensitiveDetector.h.
Referenced by createHit(), and sendHit().
float TkAccumulatingSensitiveDetector::py [private] |
Definition at line 79 of file TkAccumulatingSensitiveDetector.h.
Referenced by createHit(), and sendHit().
float TkAccumulatingSensitiveDetector::pz [private] |
Definition at line 79 of file TkAccumulatingSensitiveDetector.h.
Referenced by createHit(), and sendHit().
float TkAccumulatingSensitiveDetector::rTracker [private] |
Definition at line 90 of file TkAccumulatingSensitiveDetector.h.
Referenced by update().
Definition at line 65 of file TkAccumulatingSensitiveDetector.h.
Referenced by clearHits(), fillHits(), getNames(), sendHit(), TkAccumulatingSensitiveDetector(), and ~TkAccumulatingSensitiveDetector().
Definition at line 64 of file TkAccumulatingSensitiveDetector.h.
Referenced by clearHits(), fillHits(), getNames(), sendHit(), TkAccumulatingSensitiveDetector(), and ~TkAccumulatingSensitiveDetector().
Definition at line 73 of file TkAccumulatingSensitiveDetector.h.
Referenced by createHit(), TkAccumulatingSensitiveDetector(), and ~TkAccumulatingSensitiveDetector().
const SimTrackManager* TkAccumulatingSensitiveDetector::theManager [private] |
Definition at line 71 of file TkAccumulatingSensitiveDetector.h.
double TkAccumulatingSensitiveDetector::theSigma [private] |
Definition at line 74 of file TkAccumulatingSensitiveDetector.h.
Referenced by TkAccumulatingSensitiveDetector(), and tofBin().
float TkAccumulatingSensitiveDetector::zTracker [private] |
Definition at line 91 of file TkAccumulatingSensitiveDetector.h.
Referenced by update().