#include <MuonSensitiveDetector.h>
implementation of SensitiveDetector for the muon detector; a MuonSlaveSD handles the interfacing to the database; numbering scheme are booked according to the detector name
Modification: 19/05/03. P.Arce Add SimTracks selection
Definition at line 44 of file MuonSensitiveDetector.h.
MuonSensitiveDetector::MuonSensitiveDetector | ( | std::string | name, |
const DDCompactView & | cpv, | ||
SensitiveDetectorCatalog & | clg, | ||
edm::ParameterSet const & | p, | ||
const SimTrackManager * | manager | ||
) |
Definition at line 29 of file MuonSensitiveDetector.cc.
References SensitiveDetector::AssignSD(), detector, g4numbering, edm::ParameterSet::getParameter(), MuonSubDetector::isEndcap(), MuonSubDetector::isRpc(), LogDebug, SensitiveDetectorCatalog::logicalNames(), myG4TrackToParticleID, numbering, printHits, SensitiveDetector::Register(), slaveMuon, STallMuonsPersistent, STenergyPersistentCut, theG4ProcessTypeEnumerator, theManager, thePrinter, and theRotation.
: SensitiveTkDetector(name, cpv, clg, p), thePV(0), theHit(0), theDetUnitId(0), theTrackID(0), theManager(manager) { edm::ParameterSet m_MuonSD = p.getParameter<edm::ParameterSet>("MuonSD"); STenergyPersistentCut = m_MuonSD.getParameter<double>("EnergyThresholdForPersistency");//Default 1. GeV STallMuonsPersistent = m_MuonSD.getParameter<bool>("AllMuonsPersistent"); printHits = m_MuonSD.getParameter<bool>("PrintHits"); // // Here simply create 1 MuonSlaveSD for the moment // LogDebug("MuonSimDebug") << "create MuonSubDetector "<<name<<std::endl; detector = new MuonSubDetector(name); LogDebug("MuonSimDebug") << "create MuonFrameRotation"<<std::endl; if (detector->isEndcap()) { // cout << "MuonFrameRotation create MuonEndcapFrameRotation"<<endl; theRotation=new MuonEndcapFrameRotation(); } else if (detector->isRpc()) { // cout << "MuonFrameRotation create MuonRpcFrameRotation"<<endl; theRotation=new MuonRpcFrameRotation( cpv ); } else { theRotation = 0; } LogDebug("MuonSimDebug") << "create MuonSlaveSD"<<std::endl; slaveMuon = new MuonSlaveSD(detector,theManager); LogDebug("MuonSimDebug") << "create MuonSimHitNumberingScheme"<<std::endl; numbering = new MuonSimHitNumberingScheme(detector, cpv); g4numbering = new MuonG4Numbering(cpv); // // 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++){ LogDebug("MuonSimDebug") << name << " MuonSensitiveDetector:: attaching SD to LV " << *it << std::endl; this->AssignSD(*it); } if (printHits) { thePrinter = new SimHitPrinter("HitPositionOSCAR.dat"); } LogDebug("MuonSimDebug") << " EnergyThresholdForPersistency " << STenergyPersistentCut << " AllMuonsPersistent " << STallMuonsPersistent << std::endl; theG4ProcessTypeEnumerator = new G4ProcessTypeEnumerator; myG4TrackToParticleID = new G4TrackToParticleID; }
MuonSensitiveDetector::~MuonSensitiveDetector | ( | ) | [virtual] |
Definition at line 91 of file MuonSensitiveDetector.cc.
References detector, g4numbering, myG4TrackToParticleID, numbering, slaveMuon, theG4ProcessTypeEnumerator, and theRotation.
{ delete g4numbering; delete numbering; delete slaveMuon; delete theRotation; delete detector; delete theG4ProcessTypeEnumerator; delete myG4TrackToParticleID; }
void MuonSensitiveDetector::clearHits | ( | ) | [private, virtual] |
Implements SensitiveDetector.
Definition at line 119 of file MuonSensitiveDetector.cc.
References TrackingSlaveSD::Initialize(), LogDebug, and slaveMuon.
Referenced by update().
{ LogDebug("MuonSimDebug") << "MuonSensitiveDetector::clearHits"<<std::endl; slaveMuon->Initialize(); }
void MuonSensitiveDetector::createHit | ( | G4Step * | aStep | ) | [private] |
Definition at line 194 of file MuonSensitiveDetector.cc.
References abs, SensitiveDetector::ConvertToLocal3DPoint(), detector, SensitiveDetector::FinalStepPosition(), FinalStepPositionVsParent(), getOrCreateTrackInformation(), info, SensitiveDetector::InitialStepPosition(), InitialStepPositionVsParent(), MuonSubDetector::isBarrel(), MuonSubDetector::isEndcap(), SensitiveDetector::LocalCoordinates, LogDebug, mag(), myG4TrackToParticleID, AlCaHLTBitMon_ParallelJobs::p, p2, G4TrackToParticleID::particleID(), PV3DBase< T, PVType, FrameType >::phi(), printHits, G4ProcessTypeEnumerator::processId(), setDetUnitId(), STallMuonsPersistent, STenergyPersistentCut, TrackInformation::storeTrack(), storeVolumeAndTrack(), theDetUnitId, theG4ProcessTypeEnumerator, theGlobalEntry, theHit, thePV, PV3DBase< T, PVType, FrameType >::theta(), theTrackID, toOrcaRef(), toOrcaUnits(), csvLumiCalc::unit, SensitiveDetector::WorldCoordinates, PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().
Referenced by ProcessHits().
{ G4Track * theTrack = aStep->GetTrack(); Local3DPoint theEntryPoint; Local3DPoint theExitPoint; if (detector->isBarrel()) { theEntryPoint= toOrcaUnits(toOrcaRef(InitialStepPositionVsParent(aStep,1),aStep)); // 1 level up theExitPoint= toOrcaUnits(toOrcaRef(FinalStepPositionVsParent(aStep,1),aStep)); } else if (detector->isEndcap()) { // save local z at current level theEntryPoint= toOrcaUnits(toOrcaRef(InitialStepPosition(aStep,LocalCoordinates),aStep)); theExitPoint= toOrcaUnits(toOrcaRef(FinalStepPosition(aStep,LocalCoordinates),aStep)); float zentry = theEntryPoint.z(); float zexit = theExitPoint.z(); Local3DPoint tempEntryPoint= toOrcaUnits(toOrcaRef(InitialStepPositionVsParent(aStep,4),aStep)); // 4 levels up Local3DPoint tempExitPoint= toOrcaUnits(toOrcaRef(FinalStepPositionVsParent(aStep,4),aStep)); // reset local z from z wrt deep-parent volume to z wrt low-level volume theEntryPoint = Local3DPoint( tempEntryPoint.x(), tempEntryPoint.y(), zentry ); theExitPoint = Local3DPoint( tempExitPoint.x(), tempExitPoint.y(), zexit ); } else { theEntryPoint= toOrcaUnits(toOrcaRef(InitialStepPosition(aStep,LocalCoordinates),aStep)); theExitPoint= toOrcaUnits(toOrcaRef(FinalStepPosition(aStep,LocalCoordinates),aStep)); } float thePabs = aStep->GetPreStepPoint()->GetMomentum().mag()/GeV; float theTof = aStep->GetPreStepPoint()->GetGlobalTime()/nanosecond; float theEnergyLoss = aStep->GetTotalEnergyDeposit()/GeV; // int theParticleType = theTrack->GetDefinition()->GetPDGEncoding(); int theParticleType = myG4TrackToParticleID->particleID(theTrack); G4ThreeVector gmd = aStep->GetPreStepPoint()->GetMomentumDirection(); G4ThreeVector lmd = ((G4TouchableHistory *)(aStep->GetPreStepPoint()->GetTouchable()))->GetHistory() ->GetTopTransform().TransformAxis(gmd); Local3DPoint lnmd = toOrcaRef(ConvertToLocal3DPoint(lmd),aStep); float theThetaAtEntry = lnmd.theta(); float thePhiAtEntry = lnmd.phi(); storeVolumeAndTrack( aStep ); theDetUnitId = setDetUnitId(aStep); Global3DPoint theGlobalPos; if (printHits) { Local3DPoint theGlobalHelp = InitialStepPosition(aStep,WorldCoordinates); theGlobalEntry = toOrcaUnits(Global3DPoint (theGlobalHelp.x(),theGlobalHelp.y(),theGlobalHelp.z())); G4StepPoint * preStepPoint = aStep->GetPreStepPoint(); G4TouchableHistory * theTouchable=(G4TouchableHistory *) (preStepPoint->GetTouchable()); theGlobalHelp=ConvertToLocal3DPoint(theTouchable->GetTranslation()); theGlobalPos = toOrcaUnits(Global3DPoint (theGlobalHelp.x(),theGlobalHelp.y(),theGlobalHelp.z())); // const G4RotationMatrix * theGlobalRot = theTouchable->GetRotation(); } LogDebug("MuonSimDebug") << "MuonSensitiveDetector::createHit UpdatablePSimHit"<<std::endl; theHit = new UpdatablePSimHit(theEntryPoint,theExitPoint,thePabs,theTof, theEnergyLoss,theParticleType,theDetUnitId, theTrackID,theThetaAtEntry,thePhiAtEntry, theG4ProcessTypeEnumerator->processId(theTrack->GetCreatorProcess())); LogDebug("MuonSimDebug") <<"=== NEW ==================> ELOSS = "<<theEnergyLoss<<" " <<thePV->GetLogicalVolume()->GetName()<<std::endl; const G4VProcess* p = aStep->GetPostStepPoint()->GetProcessDefinedStep(); const G4VProcess* p2 = aStep->GetPreStepPoint()->GetProcessDefinedStep(); if (p) LogDebug("MuonSimDebug") <<" POST PROCESS = "<<p->GetProcessName()<<std::endl; if (p2) LogDebug("MuonSimDebug") <<" PRE PROCESS = "<<p2->GetProcessName()<<std::endl; LogDebug("MuonSimDebug") << "newhit theta " << theThetaAtEntry<<std::endl; LogDebug("MuonSimDebug") << "newhit phi " << thePhiAtEntry<<std::endl; LogDebug("MuonSimDebug") << "newhit pabs " << thePabs<<std::endl; LogDebug("MuonSimDebug") << "newhit tof " << theTof<<std::endl; LogDebug("MuonSimDebug") << "newhit track " << theTrackID<<std::endl; LogDebug("MuonSimDebug") << "newhit entry " << theEntryPoint<<std::endl; LogDebug("MuonSimDebug") << "newhit exit " << theExitPoint<<std::endl; LogDebug("MuonSimDebug") << "newhit eloss " << theEnergyLoss << std::endl; LogDebug("MuonSimDebug") << "newhit detid " << theDetUnitId<<std::endl; LogDebug("MuonSimDebug") << "newhit delta " << (theExitPoint-theEntryPoint)<<std::endl; LogDebug("MuonSimDebug") << "newhit deltu " << (theExitPoint-theEntryPoint).unit(); LogDebug("MuonSimDebug") << " " << (theExitPoint-theEntryPoint).mag()<<std::endl; LogDebug("MuonSimDebug") << "newhit glob " << theGlobalEntry<<std::endl; LogDebug("MuonSimDebug") << "newhit dpos " << theGlobalPos<<std::endl; LogDebug("MuonSimDebug") << "newhit drot " << std::endl; // theGlobalRot->print(LogDebug("MuonSimDebug")); // //----- SimTracks: Make it persistent? // int thePID = theTrack->GetDefinition()->GetPDGEncoding(); LogDebug("MuonSimDebug") << " checking simtrack " << thePID << " " << thePabs << " STenergyPersistentCut " << STenergyPersistentCut << std::endl; if( thePabs*GeV > STenergyPersistentCut || ( abs(thePID) == 13 && STallMuonsPersistent ) ){ TrackInformation* info = getOrCreateTrackInformation(theTrack); LogDebug("MuonSimDebug") <<" track leaving hit in muons made selected for persistency"<<std::endl; info->storeTrack(true); } }
void MuonSensitiveDetector::EndOfEvent | ( | G4HCofThisEvent * | ) | [virtual] |
Reimplemented from SensitiveDetector.
Definition at line 380 of file MuonSensitiveDetector.cc.
References saveHit().
{ // TimeMe t("MuonSensitiveDetector::EndOfEvent", false); // LogDebug("MuonSimDebug") << "MuonSensitiveDetector::EndOfEvent saving last hit en event " << std::endl; saveHit(); }
void MuonSensitiveDetector::fillHits | ( | edm::PSimHitContainer & | c, |
std::string | use | ||
) | [virtual] |
Implements SensitiveTkDetector.
Definition at line 388 of file MuonSensitiveDetector.cc.
References TrackingSlaveSD::hits(), n, TrackingSlaveSD::name(), and slaveMuon.
Local3DPoint MuonSensitiveDetector::FinalStepPositionVsParent | ( | G4Step * | currentStep, |
G4int | levelsUp | ||
) | [private] |
Definition at line 418 of file MuonSensitiveDetector.cc.
References SensitiveDetector::ConvertToLocal3DPoint().
Referenced by createHit(), and updateHit().
{ G4StepPoint * postStepPoint = currentStep->GetPostStepPoint(); G4StepPoint * preStepPoint = currentStep->GetPreStepPoint(); G4ThreeVector globalCoordinates = postStepPoint->GetPosition(); G4TouchableHistory * theTouchable = (G4TouchableHistory *) (preStepPoint->GetTouchable()); G4int depth = theTouchable->GetHistory()->GetDepth(); G4ThreeVector localCoordinates = theTouchable->GetHistory() ->GetTransform(depth-levelsUp).TransformPoint(globalCoordinates); return ConvertToLocal3DPoint(localCoordinates); }
std::vector< std::string > MuonSensitiveDetector::getNames | ( | ) | [virtual] |
Reimplemented from SensitiveDetector.
Definition at line 397 of file MuonSensitiveDetector.cc.
References TrackingSlaveSD::name(), slaveMuon, and groupFilesInBlocks::temp.
TrackInformation * MuonSensitiveDetector::getOrCreateTrackInformation | ( | const G4Track * | theTrack | ) | [private] |
Definition at line 364 of file MuonSensitiveDetector.cc.
References dtNoiseDBValidation_cfg::cerr, info, and groupFilesInBlocks::temp.
Referenced by createHit().
{ G4VUserTrackInformation* temp = gTrack->GetUserInformation(); if (temp == 0){ std::cerr <<" ERROR: no G4VUserTrackInformation available"<<std::endl; abort(); }else{ TrackInformation* info = dynamic_cast<TrackInformation*>(temp); if (info ==0){ std::cerr <<" ERROR: TkSimTrackSelection: the UserInformation does not appear to be a TrackInformation"<<std::endl; abort(); } return info; } }
const MuonSlaveSD* MuonSensitiveDetector::GetSlaveMuon | ( | ) | const [inline] |
Local3DPoint MuonSensitiveDetector::InitialStepPositionVsParent | ( | G4Step * | currentStep, |
G4int | levelsUp | ||
) | [private] |
Transform from local coordinates of a volume to local coordinates of a parent volume one or more levels up the volume hierarchy: e.g. levelsUp = 1 for immediate parent.
This is done by moving from local_1 -> global -> local_2.
Definition at line 403 of file MuonSensitiveDetector.cc.
References SensitiveDetector::ConvertToLocal3DPoint().
Referenced by createHit().
{ G4StepPoint * preStepPoint = currentStep->GetPreStepPoint(); G4ThreeVector globalCoordinates = preStepPoint->GetPosition(); G4TouchableHistory * theTouchable=(G4TouchableHistory *) (preStepPoint->GetTouchable()); G4int depth = theTouchable->GetHistory()->GetDepth(); G4ThreeVector localCoordinates = theTouchable->GetHistory() ->GetTransform(depth-levelsUp).TransformPoint(globalCoordinates); return ConvertToLocal3DPoint(localCoordinates); }
bool MuonSensitiveDetector::newHit | ( | G4Step * | aStep | ) | [private] |
Definition at line 181 of file MuonSensitiveDetector.cc.
References setDetUnitId(), lumiQTWidget::t, theDetUnitId, thePV, and theTrackID.
Referenced by ProcessHits().
{ G4VPhysicalVolume* pv = aStep->GetPreStepPoint()->GetPhysicalVolume(); G4Track * t = aStep->GetTrack(); uint32_t currentUnitId=setDetUnitId(aStep); unsigned int currentTrackID=t->GetTrackID(); //unsigned int currentEventID=G4EventManager::GetEventManager()->GetConstCurrentEvent()->GetEventID(); bool changed=((pv!=thePV) || (currentUnitId!=theDetUnitId) || (currentTrackID!=theTrackID)); return changed; }
bool MuonSensitiveDetector::ProcessHits | ( | G4Step * | aStep, |
G4TouchableHistory * | ROhist | ||
) | [virtual] |
Implements SensitiveDetector.
Definition at line 125 of file MuonSensitiveDetector.cc.
References createHit(), SensitiveDetector::InitialStepPosition(), LogDebug, newHit(), saveHit(), storeVolumeAndTrack(), updateHit(), and SensitiveDetector::WorldCoordinates.
{ LogDebug("MuonSimDebug") <<" MuonSensitiveDetector::ProcessHits "<<InitialStepPosition(aStep,WorldCoordinates)<<std::endl; // TimeMe t1( theHitTimer, false); if (aStep->GetTotalEnergyDeposit()>0.){ // do not count neutrals that are killed by User Limits MinEKine if( aStep->GetTrack()->GetDynamicParticle()->GetCharge() != 0 ){ if (newHit(aStep)) { saveHit(); createHit(aStep); } else { updateHit(aStep); } return true; } else { storeVolumeAndTrack(aStep); return false; } } return false; }
void MuonSensitiveDetector::saveHit | ( | ) | [private] |
Definition at line 343 of file MuonSensitiveDetector.cc.
References detector, PSimHit::detUnitId(), PSimHit::entryPoint(), PSimHit::exitPoint(), MuonSubDetector::name(), SimHitPrinter::printGlobal(), printHits, SimHitPrinter::printId(), SimHitPrinter::printLocal(), TrackingSlaveSD::processHits(), slaveMuon, SimHitPrinter::startNewSimHit(), theGlobalEntry, theHit, and thePrinter.
Referenced by EndOfEvent(), and ProcessHits().
{ if (theHit) { if (printHits) { thePrinter->startNewSimHit(detector->name()); thePrinter->printId(theHit->detUnitId()); // thePrinter->printTrack(theHit->trackId()); //thePrinter->printPabs(theHit->pabs()); //thePrinter->printEloss(theHit->energyLoss()); thePrinter->printLocal(theHit->entryPoint(),theHit->exitPoint()); thePrinter->printGlobal(theGlobalEntry); } slaveMuon->processHits(*theHit); // seems the hit does not want to be deleted // done by the hit collection? delete theHit; theHit = 0; //set it to 0, because you are checking that is 0 } }
uint32_t MuonSensitiveDetector::setDetUnitId | ( | G4Step * | aStep | ) | [virtual] |
Implements SensitiveDetector.
Definition at line 150 of file MuonSensitiveDetector.cc.
References MuonSimHitNumberingScheme::baseNumberToUnitNumber(), g4numbering, numbering, and MuonG4Numbering::PhysicalVolumeToBaseNumber().
Referenced by createHit(), and newHit().
{ // G4VPhysicalVolume * pv = aStep->GetPreStepPoint()->GetPhysicalVolume(); MuonBaseNumber num = g4numbering->PhysicalVolumeToBaseNumber(aStep); return numbering->baseNumberToUnitNumber(num); }
void MuonSensitiveDetector::storeVolumeAndTrack | ( | G4Step * | aStep | ) | [private] |
Definition at line 173 of file MuonSensitiveDetector.cc.
References lumiQTWidget::t, thePV, and theTrackID.
Referenced by createHit(), and ProcessHits().
{ G4VPhysicalVolume* pv = aStep->GetPreStepPoint()->GetPhysicalVolume(); G4Track * t = aStep->GetTrack(); thePV=pv; theTrackID=t->GetTrackID(); }
Local3DPoint MuonSensitiveDetector::toOrcaRef | ( | Local3DPoint | in, |
G4Step * | s | ||
) | [private] |
Definition at line 158 of file MuonSensitiveDetector.cc.
References theRotation, and MuonFrameRotation::transformPoint().
Referenced by createHit(), and updateHit().
{ if (theRotation !=0 ) { return theRotation->transformPoint(in,s); } return (in); }
Global3DPoint MuonSensitiveDetector::toOrcaUnits | ( | Global3DPoint | in | ) | [private] |
Definition at line 169 of file MuonSensitiveDetector.cc.
References PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().
{ return Global3DPoint(in.x()/cm,in.y()/cm,in.z()/cm); }
Local3DPoint MuonSensitiveDetector::toOrcaUnits | ( | Local3DPoint | in | ) | [private] |
Definition at line 165 of file MuonSensitiveDetector.cc.
References PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().
Referenced by createHit(), and updateHit().
{ return Local3DPoint(in.x()/cm,in.y()/cm,in.z()/cm); }
std::string MuonSensitiveDetector::type | ( | ) |
void MuonSensitiveDetector::update | ( | const ::EndOfEvent * | ev | ) | [private] |
Definition at line 113 of file MuonSensitiveDetector.cc.
{
//slaveMuon->renumbering(theManager);
}
void MuonSensitiveDetector::update | ( | const BeginOfEvent * | ) | [private, virtual] |
This routine will be called when the appropriate signal arrives.
Implements Observer< const BeginOfEvent * >.
Definition at line 103 of file MuonSensitiveDetector.cc.
References clearHits(), theDetUnitId, thePV, and theTrackID.
{ clearHits(); //----- Initialize variables to check if two steps belong to same hit thePV = 0; theDetUnitId = 0; theTrackID = 0; }
void MuonSensitiveDetector::updateHit | ( | G4Step * | aStep | ) | [private] |
Definition at line 298 of file MuonSensitiveDetector.cc.
References UpdatablePSimHit::addEnergyLoss(), dtNoiseDBValidation_cfg::cerr, detector, PSimHit::energyLoss(), PSimHit::entryPoint(), SensitiveDetector::FinalStepPosition(), FinalStepPositionVsParent(), MuonSubDetector::isBarrel(), MuonSubDetector::isEndcap(), SensitiveDetector::LocalCoordinates, LogDebug, mag(), AlCaHLTBitMon_ParallelJobs::p, p2, theDetUnitId, theHit, thePV, toOrcaRef(), toOrcaUnits(), csvLumiCalc::unit, UpdatablePSimHit::updateExitPoint(), PV3DBase< T, PVType, FrameType >::x(), PV3DBase< T, PVType, FrameType >::y(), and PV3DBase< T, PVType, FrameType >::z().
Referenced by ProcessHits().
{ // float thePabs = aStep->GetPreStepPoint()->GetMomentum().mag()/GeV; // Local3DPoint theEntryPoint= InitialStepPosition(aStep,LocalCoordinates); Local3DPoint theExitPoint; if (detector->isBarrel()) { theExitPoint= toOrcaUnits(toOrcaRef(FinalStepPositionVsParent(aStep,1),aStep)); } else if (detector->isEndcap()) { // save local z at current level theExitPoint= toOrcaUnits(toOrcaRef(FinalStepPosition(aStep,LocalCoordinates),aStep)); float zexit = theExitPoint.z(); Local3DPoint tempExitPoint= toOrcaUnits(toOrcaRef(FinalStepPositionVsParent(aStep,4),aStep)); theExitPoint = Local3DPoint( tempExitPoint.x(), tempExitPoint.y(), zexit ); } else { theExitPoint= toOrcaUnits(toOrcaRef(FinalStepPosition(aStep,LocalCoordinates),aStep)); } float theEnergyLoss = aStep->GetTotalEnergyDeposit()/GeV; if( theHit == 0 ){ std::cerr << "!!ERRROR in MuonSensitiveDetector::updateHit. It is called when there is no hit " << std::endl; } theHit->updateExitPoint(theExitPoint); theHit->addEnergyLoss(theEnergyLoss); LogDebug("MuonSimDebug") <<"=== UPDATE ===============> ELOSS = "<<theEnergyLoss<<" " <<thePV->GetLogicalVolume()->GetName()<<std::endl; const G4VProcess* p = aStep->GetPostStepPoint()->GetProcessDefinedStep(); const G4VProcess* p2 = aStep->GetPreStepPoint()->GetProcessDefinedStep(); if (p) LogDebug("MuonSimDebug") <<" POST PROCESS = "<<p->GetProcessName()<<std::endl; if (p2) LogDebug("MuonSimDebug") <<" PRE PROCESS = "<<p2->GetProcessName()<<std::endl; LogDebug("MuonSimDebug") << "updhit exit " << theExitPoint<<std::endl; LogDebug("MuonSimDebug") << "updhit eloss " << theHit->energyLoss() <<std::endl; LogDebug("MuonSimDebug") << "updhit detid " << theDetUnitId<<std::endl; LogDebug("MuonSimDebug") << "updhit delta " << (theExitPoint-theHit->entryPoint())<<std::endl; LogDebug("MuonSimDebug") << "updhit deltu " << (theExitPoint-theHit->entryPoint()).unit(); LogDebug("MuonSimDebug") << " " << (theExitPoint-theHit->entryPoint()).mag()<<std::endl; }
MuonSubDetector* MuonSensitiveDetector::detector [private] |
Definition at line 80 of file MuonSensitiveDetector.h.
Referenced by createHit(), MuonSensitiveDetector(), saveHit(), updateHit(), and ~MuonSensitiveDetector().
Definition at line 82 of file MuonSensitiveDetector.h.
Referenced by MuonSensitiveDetector(), setDetUnitId(), and ~MuonSensitiveDetector().
Definition at line 113 of file MuonSensitiveDetector.h.
Referenced by createHit(), MuonSensitiveDetector(), and ~MuonSensitiveDetector().
Definition at line 79 of file MuonSensitiveDetector.h.
Referenced by MuonSensitiveDetector(), setDetUnitId(), and ~MuonSensitiveDetector().
bool MuonSensitiveDetector::printHits [private] |
Definition at line 103 of file MuonSensitiveDetector.h.
Referenced by createHit(), MuonSensitiveDetector(), and saveHit().
MuonSlaveSD* MuonSensitiveDetector::slaveMuon [private] |
Definition at line 78 of file MuonSensitiveDetector.h.
Referenced by clearHits(), fillHits(), getNames(), GetSlaveMuon(), MuonSensitiveDetector(), saveHit(), and ~MuonSensitiveDetector().
bool MuonSensitiveDetector::STallMuonsPersistent [private] |
Definition at line 109 of file MuonSensitiveDetector.h.
Referenced by createHit(), and MuonSensitiveDetector().
double MuonSensitiveDetector::STenergyPersistentCut [private] |
Definition at line 108 of file MuonSensitiveDetector.h.
Referenced by createHit(), and MuonSensitiveDetector().
uint32_t MuonSensitiveDetector::theDetUnitId [private] |
Definition at line 100 of file MuonSensitiveDetector.h.
Referenced by createHit(), newHit(), update(), and updateHit().
Definition at line 111 of file MuonSensitiveDetector.h.
Referenced by createHit(), MuonSensitiveDetector(), and ~MuonSensitiveDetector().
Definition at line 105 of file MuonSensitiveDetector.h.
Referenced by createHit(), and saveHit().
UpdatablePSimHit* MuonSensitiveDetector::theHit [private] |
Definition at line 99 of file MuonSensitiveDetector.h.
Referenced by createHit(), saveHit(), and updateHit().
const SimTrackManager* MuonSensitiveDetector::theManager [private] |
Definition at line 114 of file MuonSensitiveDetector.h.
Referenced by MuonSensitiveDetector().
SimHitPrinter* MuonSensitiveDetector::thePrinter [private] |
Definition at line 104 of file MuonSensitiveDetector.h.
Referenced by MuonSensitiveDetector(), and saveHit().
G4VPhysicalVolume* MuonSensitiveDetector::thePV [private] |
Definition at line 98 of file MuonSensitiveDetector.h.
Referenced by createHit(), newHit(), storeVolumeAndTrack(), update(), and updateHit().
Definition at line 81 of file MuonSensitiveDetector.h.
Referenced by MuonSensitiveDetector(), toOrcaRef(), and ~MuonSensitiveDetector().
unsigned int MuonSensitiveDetector::theTrackID [private] |
Definition at line 101 of file MuonSensitiveDetector.h.
Referenced by createHit(), newHit(), storeVolumeAndTrack(), and update().