#include <FP420SD.h>
Public Member Functions | |
virtual void | clear () |
virtual void | DrawAll () |
virtual void | EndOfEvent (G4HCofThisEvent *eventHC) |
void | fillHits (edm::PSimHitContainer &, std::string use) |
FP420SD (std::string, const DDCompactView &, SensitiveDetectorCatalog &, edm::ParameterSet const &, const SimTrackManager *) | |
virtual double | getEnergyDeposit (G4Step *step) |
std::vector< std::string > | getNames () |
virtual void | Initialize (G4HCofThisEvent *HCE) |
virtual void | PrintAll () |
virtual bool | ProcessHits (G4Step *, G4TouchableHistory *) |
virtual uint32_t | setDetUnitId (G4Step *) |
virtual | ~FP420SD () |
Protected Attributes | |
float | edepositEM |
float | edepositHAD |
G4int | emPDG |
G4int | epPDG |
G4int | gammaPDG |
Private Member Functions | |
virtual void | clearHits () |
void | CreateNewHit () |
void | GetStepInfo (G4Step *aStep) |
G4bool | HitExists () |
void | ResetForNewPrimary () |
G4ThreeVector | SetToLocal (G4ThreeVector global) |
G4ThreeVector | SetToLocalExit (G4ThreeVector globalPoint) |
void | StoreHit (FP420G4Hit *) |
void | Summarize () |
void | update (const BeginOfRun *) |
This routine will be called when the appropriate signal arrives. | |
void | update (const ::EndOfEvent *) |
void | update (const BeginOfEvent *) |
This routine will be called when the appropriate signal arrives. | |
void | UpdateHit () |
Private Attributes | |
FP420G4Hit * | currentHit |
G4VPhysicalVolume * | currentPV |
float | edeposit |
float | Eloss |
G4ThreeVector | entrancePoint |
int | eventno |
G4ThreeVector | exitPoint |
G4int | hcID |
G4ThreeVector | hitPoint |
G4ThreeVector | hitPointExit |
G4ThreeVector | hitPointLocal |
G4ThreeVector | hitPointLocalExit |
float | incidentEnergy |
std::string | name |
FP420NumberingScheme * | numberingScheme |
float | Pabs |
int | ParentId |
short | ParticleType |
float | PhiAtEntry |
G4StepPoint * | postStepPoint |
G4StepPoint * | preStepPoint |
uint32_t | previousUnitID |
unsigned int | primaryID |
unsigned int | primID |
TrackingSlaveSD * | slave |
G4ThreeVector | theEntryPoint |
G4ThreeVector | theExitPoint |
FP420G4HitCollection * | theHC |
const SimTrackManager * | theManager |
float | ThetaAtEntry |
G4Track * | theTrack |
float | Tof |
G4int | tsID |
G4double | tSlice |
G4int | tSliceID |
uint32_t | unitID |
float | Vx |
float | Vy |
float | Vz |
float | X |
float | Y |
float | Z |
FP420SD::FP420SD | ( | std::string | name, |
const DDCompactView & | cpv, | ||
SensitiveDetectorCatalog & | clg, | ||
edm::ParameterSet const & | p, | ||
const SimTrackManager * | manager | ||
) |
Definition at line 55 of file FP420SD.cc.
References SensitiveDetector::AssignSD(), gather_cfg::cout, edm::ParameterSet::getParameter(), edm::ParameterSet::getUntrackedParameter(), LogDebug, SensitiveDetectorCatalog::logicalNames(), numberingScheme, SensitiveDetector::Register(), and slave.
: SensitiveTkDetector(name, cpv, clg, p), numberingScheme(0), name(name), hcID(-1), theHC(0), theManager(manager), currentHit(0), theTrack(0), currentPV(0), unitID(0), previousUnitID(0), preStepPoint(0), postStepPoint(0), eventno(0){ //------------------------------------------------------------------- /* FP420SD::FP420SD(G4String name, const DDCompactView & cpv, edm::ParameterSet const & p, const SimTrackManager* manager) : CaloSD(name, cpv, p, manager), numberingScheme(0), name(name),hcID(-1), theHC(0), currentHit(0), theTrack(0), currentPV(0), unitID(0), previousUnitID(0), preStepPoint(0), postStepPoint(0), eventno(0){ */ //------------------------------------------------------------------- //Add FP420 Sentitive Detector Name collectionName.insert(name); //Parameters edm::ParameterSet m_p = p.getParameter<edm::ParameterSet>("FP420SD"); int verbn = m_p.getUntrackedParameter<int>("Verbosity"); //int verbn = 1; SetVerboseLevel(verbn); LogDebug("FP420Sim") << "*******************************************************\n" << "* *\n" << "* Constructing a FP420SD with name " << name << "\n" << "* *\n" << "*******************************************************"; slave = new TrackingSlaveSD(name); // // attach detectors (LogicalVolumes) // std::vector<std::string> lvNames = clg.logicalNames(name); this->Register(); for (std::vector<std::string>::iterator it=lvNames.begin(); it !=lvNames.end(); it++) { this->AssignSD(*it); edm::LogInfo("FP420Sim") << "FP420SD : Assigns SD to LV " << (*it); } if (name == "FP420SI") { if (verbn > 0) { edm::LogInfo("FP420Sim") << "name = FP420SI and new FP420NumberingSchem"; std::cout << "name = FP420SI and new FP420NumberingSchem"<< std::endl; } numberingScheme = new FP420NumberingScheme() ; } else { edm::LogWarning("FP420Sim") << "FP420SD: ReadoutName not supported\n"; std::cout << "FP420SD: ReadoutName not supported"<< std::endl; } edm::LogInfo("FP420Sim") << "FP420SD: Instantiation completed"; std::cout << "FP420SD: Instantiation completed"<< std::endl; }
FP420SD::~FP420SD | ( | ) | [virtual] |
Definition at line 125 of file FP420SD.cc.
References numberingScheme, and slave.
{ //AZ: if (slave) delete slave; if (numberingScheme) delete numberingScheme; }
void FP420SD::clear | ( | void | ) | [virtual] |
Definition at line 485 of file FP420SD.cc.
{ }
void FP420SD::clearHits | ( | ) | [private, virtual] |
Implements SensitiveDetector.
Definition at line 534 of file FP420SD.cc.
References TrackingSlaveSD::Initialize(), and slave.
Referenced by update().
{ //AZ: slave->Initialize(); }
void FP420SD::CreateNewHit | ( | ) | [private] |
Definition at line 297 of file FP420SD.cc.
References currentHit, currentPV, Eloss, hitPoint, hitPointLocal, hitPointLocalExit, incidentEnergy, LogDebug, NULL, Pabs, ParentId, ParticleType, PhiAtEntry, previousUnitID, primaryID, primID, FP420G4Hit::setEnergyLoss(), FP420G4Hit::setEntry(), FP420G4Hit::setEntryLocalP(), FP420G4Hit::setExitLocalP(), FP420G4Hit::setIncidentEnergy(), FP420G4Hit::setPabs(), FP420G4Hit::setParentId(), FP420G4Hit::setParticleType(), FP420G4Hit::setPhiAtEntry(), FP420G4Hit::setThetaAtEntry(), FP420G4Hit::setTimeSlice(), FP420G4Hit::setTof(), FP420G4Hit::setTrackID(), FP420G4Hit::setUnitID(), FP420G4Hit::setVx(), FP420G4Hit::setVy(), FP420G4Hit::setVz(), FP420G4Hit::setX(), FP420G4Hit::setY(), FP420G4Hit::setZ(), StoreHit(), ThetaAtEntry, theTrack, Tof, tsID, tSlice, tSliceID, unitID, Vx, Vy, Vz, X, Y, and Z.
Referenced by ProcessHits().
{ #ifdef debug // << " MVid = " << currentPV->GetMother()->GetCopyNo() LogDebug("FP420Sim") << "FP420SD CreateNewHit for" << " PV " << currentPV->GetName() << " PVid = " << currentPV->GetCopyNo() << " Unit " << unitID <<std::endl; LogDebug("FP420Sim") << " primary " << primaryID << " time slice " << tSliceID << " For Track " << theTrack->GetTrackID() << " which is a " << theTrack->GetDefinition()->GetParticleName(); if (theTrack->GetTrackID()==1) { LogDebug("FP420Sim") << " of energy " << theTrack->GetTotalEnergy(); } else { LogDebug("FP420Sim") << " daughter of part. " << theTrack->GetParentID(); } LogDebug("FP420Sim") << " and created by " ; if (theTrack->GetCreatorProcess()!=NULL) LogDebug("FP420Sim") << theTrack->GetCreatorProcess()->GetProcessName() ; else LogDebug("FP420Sim") << "NO process"; LogDebug("FP420Sim") << std::endl; #endif currentHit = new FP420G4Hit; currentHit->setTrackID(primaryID); currentHit->setTimeSlice(tSlice); currentHit->setUnitID(unitID); currentHit->setIncidentEnergy(incidentEnergy); currentHit->setPabs(Pabs); currentHit->setTof(Tof); currentHit->setEnergyLoss(Eloss); currentHit->setParticleType(ParticleType); currentHit->setThetaAtEntry(ThetaAtEntry); currentHit->setPhiAtEntry(PhiAtEntry); // currentHit->setEntry(entrancePoint); currentHit->setEntry(hitPoint); currentHit->setEntryLocalP(hitPointLocal); currentHit->setExitLocalP(hitPointLocalExit); currentHit->setParentId(ParentId); currentHit->setVx(Vx); currentHit->setVy(Vy); currentHit->setVz(Vz); currentHit->setX(X); currentHit->setY(Y); currentHit->setZ(Z); //AZ:12.10.2007 // UpdateHit(); // buffer for next steps: tsID = tSliceID; primID = primaryID; previousUnitID = unitID; StoreHit(currentHit); }
void FP420SD::DrawAll | ( | ) | [virtual] |
Definition at line 489 of file FP420SD.cc.
{ }
void FP420SD::EndOfEvent | ( | G4HCofThisEvent * | eventHC | ) | [virtual] |
Reimplemented from SensitiveDetector.
Definition at line 403 of file FP420SD.cc.
References FP420G4Hit::getEnergyLoss(), FP420G4Hit::getEntry(), FP420G4Hit::getEntryLocalP(), FP420G4Hit::getExitLocalP(), FP420G4Hit::getPabs(), FP420G4Hit::getParticleType(), FP420G4Hit::getPhiAtEntry(), FP420G4Hit::getThetaAtEntry(), FP420G4Hit::getTof(), FP420G4Hit::getTrackID(), FP420G4Hit::getUnitID(), j, LogDebug, TrackingSlaveSD::processHits(), slave, Summarize(), and theHC.
{ // here we loop over transient hits and make them persistent // if(theHC->entries() > 100){ // LogDebug("FP420Sim") << "FP420SD: warning!!! Number of hits exceed 100 and =" << theHC->entries() << "\n"; // } // for (int j=0; j<theHC->entries() && j<100; j++) { int nhitsHPS240=0; int nhitsFP420=0; for (int j=0; j<theHC->entries(); j++) { FP420G4Hit* aHit = (*theHC)[j]; if((fabs(aHit->getTof())> 780. && fabs(aHit->getTof())< 840.)) ++nhitsHPS240; if((fabs(aHit->getTof())>1380. && fabs(aHit->getTof())<1450.)) ++nhitsFP420; // if(fabs(aHit->getTof()) < 1700.) { if((fabs(aHit->getTof())>780. && fabs(aHit->getTof())<840. && nhitsHPS240<200.) || ( fabs(aHit->getTof())>1380. && fabs(aHit->getTof())<1450. && nhitsFP420<200.)) { #ifdef ddebug // LogDebug("FP420SD") << " FP420Hit " << j << " " << *aHit << std::endl; LogDebug("FP420Sim") << "hit number" << j << "unit ID = "<<aHit->getUnitID()<< "\n"; LogDebug("FP420Sim") << "entry z " << aHit->getEntry().z()<< "\n"; LogDebug("FP420Sim") << "entr theta " << aHit->getThetaAtEntry()<< "\n"; #endif // Local3DPoint locExitPoint(0,0,0); // Local3DPoint locEntryPoint(aHit->getEntry().x(), // aHit->getEntry().y(), // aHit->getEntry().z()); Local3DPoint locExitPoint(aHit->getExitLocalP().x(), aHit->getExitLocalP().y(), aHit->getExitLocalP().z()); Local3DPoint locEntryPoint(aHit->getEntryLocalP().x(), aHit->getEntryLocalP().y(), aHit->getEntryLocalP().z()); // implicit conversion (slicing) to PSimHit!!! slave->processHits(PSimHit(locEntryPoint,locExitPoint,//entryPoint(), exitPoint() Local3DPoint aHit->getPabs(),// pabs() float aHit->getTof(), // tof() float aHit->getEnergyLoss(), // energyLoss() float aHit->getParticleType(),// particleType() int aHit->getUnitID(), // detUnitId() unsigned int aHit->getTrackID(),// trackId() unsigned int aHit->getThetaAtEntry(),// thetaAtEntry() float aHit->getPhiAtEntry())); // phiAtEntry() float //PSimHit( const Local3DPoint& entry, const Local3DPoint& exit, // float pabs, float tof, float eloss, int particleType, // unsigned int detId, unsigned int trackId, // float theta, float phi, unsigned short processType=0) : // LocalVector direction = hit.exitPoint() - hit.entryPoint(); //hit.energyLoss /* aHit->getEM(), - aHit->getHadr(), - aHit->getIncidentEnergy(), - aHit->getTimeSlice(), - aHit->getEntry(), - aHit->getParentId(), aHit->getEntryLocalP(), - aHit->getExitLocalP(), - aHit->getX(), - aHit->getY(), - aHit->getZ(), - aHit->getVx(), - aHit->getVy(), - aHit->getVz())); - */ }//if Tof<1600. if nhits<100 } // for loop on hits Summarize(); }
void FP420SD::fillHits | ( | edm::PSimHitContainer & | c, |
std::string | use | ||
) | [virtual] |
Implements SensitiveTkDetector.
Definition at line 507 of file FP420SD.cc.
References gather_cfg::cout, TrackingSlaveSD::hits(), n, name, TrackingSlaveSD::name(), and slave.
double FP420SD::getEnergyDeposit | ( | G4Step * | step | ) | [virtual] |
Definition at line 134 of file FP420SD.cc.
Referenced by GetStepInfo().
{
return aStep->GetTotalEnergyDeposit();
}
std::vector< std::string > FP420SD::getNames | ( | ) | [virtual] |
Reimplemented from SensitiveDetector.
Definition at line 539 of file FP420SD.cc.
References TrackingSlaveSD::name(), slave, and groupFilesInBlocks::temp.
void FP420SD::GetStepInfo | ( | G4Step * | aStep | ) | [private] |
Definition at line 177 of file FP420SD.cc.
References currentPV, edeposit, edepositEM, edepositHAD, Eloss, emPDG, epPDG, gammaPDG, getEnergyDeposit(), hitPoint, hitPointExit, hitPointLocal, hitPointLocalExit, LogDebug, Pabs, ParentId, ParticleType, PhiAtEntry, postStepPoint, preStepPoint, primaryID, setDetUnitId(), ThetaAtEntry, theTrack, Tof, tSlice, tSliceID, unitID, Vx, Vy, Vz, X, Y, and Z.
Referenced by ProcessHits().
{ preStepPoint = aStep->GetPreStepPoint(); postStepPoint= aStep->GetPostStepPoint(); theTrack = aStep->GetTrack(); hitPoint = preStepPoint->GetPosition(); currentPV = preStepPoint->GetPhysicalVolume(); hitPointExit = postStepPoint->GetPosition(); hitPointLocal = preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(hitPoint); hitPointLocalExit = preStepPoint->GetTouchable()->GetHistory()->GetTopTransform().TransformPoint(hitPointExit); G4int particleCode = theTrack->GetDefinition()->GetPDGEncoding(); if (particleCode == emPDG || particleCode == epPDG || particleCode == gammaPDG ) { edepositEM = getEnergyDeposit(aStep); edepositHAD = 0.; } else { edepositEM = 0.; edepositHAD = getEnergyDeposit(aStep); } edeposit = aStep->GetTotalEnergyDeposit(); tSlice = (postStepPoint->GetGlobalTime() )/nanosecond; tSliceID = (int) tSlice; unitID = setDetUnitId(aStep); #ifdef debug LogDebug("FP420Sim") << "unitID=" << unitID <<std::endl; #endif primaryID = theTrack->GetTrackID(); // Position = hitPoint; Pabs = aStep->GetPreStepPoint()->GetMomentum().mag()/GeV; //Tof = 1400. + aStep->GetPostStepPoint()->GetGlobalTime()/nanosecond; Tof = aStep->GetPostStepPoint()->GetGlobalTime()/nanosecond; Eloss = aStep->GetTotalEnergyDeposit()/GeV; ParticleType = theTrack->GetDefinition()->GetPDGEncoding(); ThetaAtEntry = aStep->GetPreStepPoint()->GetPosition().theta()/deg; PhiAtEntry = aStep->GetPreStepPoint()->GetPosition().phi()/deg; ParentId = theTrack->GetParentID(); Vx = theTrack->GetVertexPosition().x(); Vy = theTrack->GetVertexPosition().y(); Vz = theTrack->GetVertexPosition().z(); X = hitPoint.x(); Y = hitPoint.y(); Z = hitPoint.z(); }
G4bool FP420SD::HitExists | ( | ) | [private] |
Definition at line 230 of file FP420SD.cc.
References currentHit, newFWLiteAna::found, FP420G4Hit::getTimeSliceID(), FP420G4Hit::getTrackID(), FP420G4Hit::getUnitID(), j, previousUnitID, primaryID, primID, ResetForNewPrimary(), theHC, tsID, tSliceID, unitID, and UpdateHit().
Referenced by ProcessHits().
{ if (primaryID<1) { edm::LogWarning("FP420Sim") << "***** FP420SD error: primaryID = " << primaryID << " maybe detector name changed"; } // Update if in the same detector, time-slice and for same track // if (primaryID == primID && tSliceID == tsID && unitID==previousUnitID) { if (tSliceID == tsID && unitID==previousUnitID) { //AZ: UpdateHit(); return true; } // Reset entry point for new primary if (primaryID != primID) ResetForNewPrimary(); //look in the HitContainer whether a hit with the same primID, unitID, //tSliceID already exists: G4bool found = false; // LogDebug("FP420Sim") << "FP420SD: HCollection= " << theHC->entries() <<std::endl; for (int j=0; j<theHC->entries()&&!found; j++) { FP420G4Hit* aPreviousHit = (*theHC)[j]; if (aPreviousHit->getTrackID() == primaryID && aPreviousHit->getTimeSliceID() == tSliceID && aPreviousHit->getUnitID() == unitID ) { //AZ: currentHit = aPreviousHit; found = true; } } if (found) { //AZ: UpdateHit(); return true; } else { return false; } }
void FP420SD::Initialize | ( | G4HCofThisEvent * | HCE | ) | [virtual] |
Reimplemented from SensitiveDetector.
Definition at line 138 of file FP420SD.cc.
References hcID, LogDebug, name, primID, theHC, and tsID.
{ #ifdef debug LogDebug("FP420Sim") << "FP420SD : Initialize called for " << name << std::endl; #endif theHC = new FP420G4HitCollection(name, collectionName[0]); if (hcID<0) hcID = G4SDManager::GetSDMpointer()->GetCollectionID(collectionName[0]); HCE->AddHitsCollection(hcID, theHC); tsID = -2; // primID = -2; primID = 0; }
void FP420SD::PrintAll | ( | ) | [virtual] |
bool FP420SD::ProcessHits | ( | G4Step * | aStep, |
G4TouchableHistory * | |||
) | [virtual] |
Implements SensitiveDetector.
Definition at line 156 of file FP420SD.cc.
References CreateNewHit(), edeposit, GetStepInfo(), HitExists(), LogDebug, NULL, and theHC.
{ if (aStep == NULL) { return true; } else { GetStepInfo(aStep); // LogDebug("FP420Sim") << edeposit <<std::endl; //AZ #ifdef debug LogDebug("FP420Sim") << "FP420SD : number of hits = " << theHC->entries() << std::endl; #endif if (HitExists() == false && edeposit>0. && theHC->entries()< 200 ){ CreateNewHit(); return true; } } return true; }
void FP420SD::ResetForNewPrimary | ( | ) | [private] |
Definition at line 276 of file FP420SD.cc.
References entrancePoint, exitPoint, hitPoint, hitPointExit, incidentEnergy, preStepPoint, SetToLocal(), and SetToLocalExit().
Referenced by HitExists().
{ entrancePoint = SetToLocal(hitPoint); exitPoint = SetToLocalExit(hitPointExit); incidentEnergy = preStepPoint->GetKineticEnergy(); }
uint32_t FP420SD::setDetUnitId | ( | G4Step * | aStep | ) | [virtual] |
Implements SensitiveDetector.
Definition at line 224 of file FP420SD.cc.
References FP420NumberingScheme::getUnitID(), and numberingScheme.
Referenced by GetStepInfo().
{ return (numberingScheme == 0 ? 0 : numberingScheme->getUnitID(aStep)); }
G4ThreeVector FP420SD::SetToLocal | ( | G4ThreeVector | global | ) | [private] |
Definition at line 387 of file FP420SD.cc.
References preStepPoint, and theEntryPoint.
Referenced by ResetForNewPrimary().
{ const G4VTouchable* touch= preStepPoint->GetTouchable(); theEntryPoint = touch->GetHistory()->GetTopTransform().TransformPoint(global); return theEntryPoint; }
G4ThreeVector FP420SD::SetToLocalExit | ( | G4ThreeVector | globalPoint | ) | [private] |
Definition at line 395 of file FP420SD.cc.
References postStepPoint, and theExitPoint.
Referenced by ResetForNewPrimary().
{ const G4VTouchable* touch= postStepPoint->GetTouchable(); theExitPoint = touch->GetHistory()->GetTopTransform().TransformPoint(globalPoint); return theExitPoint; }
void FP420SD::StoreHit | ( | FP420G4Hit * | hit | ) | [private] |
Definition at line 285 of file FP420SD.cc.
References theHC.
Referenced by CreateNewHit().
{ // if (primID<0) return; if (hit == 0 ) { edm::LogWarning("FP420Sim") << "FP420SD: hit to be stored is NULL !!"; return; } theHC->insert( hit ); }
void FP420SD::Summarize | ( | ) | [private] |
void FP420SD::update | ( | const ::EndOfEvent * | ) | [private] |
Definition at line 531 of file FP420SD.cc.
{ }
void FP420SD::update | ( | const BeginOfEvent * | ) | [private, virtual] |
This routine will be called when the appropriate signal arrives.
Implements Observer< const BeginOfEvent * >.
Definition at line 514 of file FP420SD.cc.
References clearHits(), eventno, and LogDebug.
void FP420SD::update | ( | const BeginOfRun * | ) | [private, virtual] |
This routine will be called when the appropriate signal arrives.
Implements Observer< const BeginOfRun * >.
Definition at line 521 of file FP420SD.cc.
References emPDG, epPDG, and gammaPDG.
{ G4ParticleTable * theParticleTable = G4ParticleTable::GetParticleTable(); G4String particleName; emPDG = theParticleTable->FindParticle(particleName="e-")->GetPDGEncoding(); epPDG = theParticleTable->FindParticle(particleName="e+")->GetPDGEncoding(); gammaPDG = theParticleTable->FindParticle(particleName="gamma")->GetPDGEncoding(); }
void FP420SD::UpdateHit | ( | ) | [private] |
Definition at line 363 of file FP420SD.cc.
References FP420G4Hit::addEnergyDeposit(), FP420G4Hit::addEnergyLoss(), currentHit, edepositEM, edepositHAD, Eloss, LogDebug, postStepPoint, previousUnitID, primaryID, primID, tsID, tSliceID, and unitID.
Referenced by HitExists().
{ // if (Eloss > 0.) { currentHit->addEnergyDeposit(edepositEM,edepositHAD); #ifdef debug LogDebug("FP420Sim") << "updateHit: add eloss " << Eloss <<std::endl; LogDebug("FP420Sim") << "CurrentHit=" << currentHit << ", PostStepPoint=" << postStepPoint->GetPosition() << std::endl; #endif //AZ // currentHit->setEnergyLoss(Eloss); currentHit->addEnergyLoss(Eloss); } // buffer for next steps: tsID = tSliceID; primID = primaryID; previousUnitID = unitID; }
FP420G4Hit* FP420SD::currentHit [private] |
Definition at line 143 of file FP420SD.h.
Referenced by CreateNewHit(), HitExists(), and UpdateHit().
G4VPhysicalVolume* FP420SD::currentPV [private] |
Definition at line 145 of file FP420SD.h.
Referenced by CreateNewHit(), and GetStepInfo().
float FP420SD::edeposit [private] |
Definition at line 155 of file FP420SD.h.
Referenced by GetStepInfo(), and ProcessHits().
float FP420SD::edepositEM [protected] |
Definition at line 182 of file FP420SD.h.
Referenced by GetStepInfo(), and UpdateHit().
float FP420SD::edepositHAD [protected] |
Definition at line 182 of file FP420SD.h.
Referenced by GetStepInfo(), and UpdateHit().
float FP420SD::Eloss [private] |
Definition at line 164 of file FP420SD.h.
Referenced by CreateNewHit(), GetStepInfo(), and UpdateHit().
G4int FP420SD::emPDG [protected] |
Definition at line 183 of file FP420SD.h.
Referenced by GetStepInfo(), and update().
G4ThreeVector FP420SD::entrancePoint [private] |
Definition at line 125 of file FP420SD.h.
Referenced by ResetForNewPrimary().
G4int FP420SD::epPDG [protected] |
Definition at line 184 of file FP420SD.h.
Referenced by GetStepInfo(), and update().
int FP420SD::eventno [private] |
G4ThreeVector FP420SD::exitPoint [private] |
Definition at line 125 of file FP420SD.h.
Referenced by ResetForNewPrimary().
G4int FP420SD::gammaPDG [protected] |
Definition at line 185 of file FP420SD.h.
Referenced by GetStepInfo(), and update().
G4int FP420SD::hcID [private] |
Definition at line 138 of file FP420SD.h.
Referenced by Initialize().
G4ThreeVector FP420SD::hitPoint [private] |
Definition at line 157 of file FP420SD.h.
Referenced by CreateNewHit(), GetStepInfo(), and ResetForNewPrimary().
G4ThreeVector FP420SD::hitPointExit [private] |
Definition at line 159 of file FP420SD.h.
Referenced by GetStepInfo(), and ResetForNewPrimary().
G4ThreeVector FP420SD::hitPointLocal [private] |
Definition at line 160 of file FP420SD.h.
Referenced by CreateNewHit(), and GetStepInfo().
G4ThreeVector FP420SD::hitPointLocalExit [private] |
Definition at line 161 of file FP420SD.h.
Referenced by CreateNewHit(), and GetStepInfo().
float FP420SD::incidentEnergy [private] |
Definition at line 134 of file FP420SD.h.
Referenced by CreateNewHit(), and ResetForNewPrimary().
std::string FP420SD::name [private] |
Reimplemented from SensitiveDetector.
Definition at line 137 of file FP420SD.h.
Referenced by fillHits(), and Initialize().
FP420NumberingScheme* FP420SD::numberingScheme [private] |
Definition at line 123 of file FP420SD.h.
Referenced by FP420SD(), setDetUnitId(), and ~FP420SD().
float FP420SD::Pabs [private] |
Definition at line 162 of file FP420SD.h.
Referenced by CreateNewHit(), and GetStepInfo().
int FP420SD::ParentId [private] |
Definition at line 170 of file FP420SD.h.
Referenced by CreateNewHit(), and GetStepInfo().
short FP420SD::ParticleType [private] |
Definition at line 165 of file FP420SD.h.
Referenced by CreateNewHit(), and GetStepInfo().
float FP420SD::PhiAtEntry [private] |
Definition at line 168 of file FP420SD.h.
Referenced by CreateNewHit(), and GetStepInfo().
G4StepPoint* FP420SD::postStepPoint [private] |
Definition at line 154 of file FP420SD.h.
Referenced by GetStepInfo(), SetToLocalExit(), and UpdateHit().
G4StepPoint* FP420SD::preStepPoint [private] |
Definition at line 153 of file FP420SD.h.
Referenced by GetStepInfo(), ResetForNewPrimary(), and SetToLocal().
uint32_t FP420SD::previousUnitID [private] |
Definition at line 147 of file FP420SD.h.
Referenced by CreateNewHit(), HitExists(), and UpdateHit().
unsigned int FP420SD::primaryID [private] |
Definition at line 149 of file FP420SD.h.
Referenced by CreateNewHit(), GetStepInfo(), HitExists(), and UpdateHit().
unsigned int FP420SD::primID [private] |
Definition at line 149 of file FP420SD.h.
Referenced by CreateNewHit(), HitExists(), Initialize(), and UpdateHit().
TrackingSlaveSD* FP420SD::slave [private] |
Definition at line 122 of file FP420SD.h.
Referenced by clearHits(), EndOfEvent(), fillHits(), FP420SD(), getNames(), and ~FP420SD().
G4ThreeVector FP420SD::theEntryPoint [private] |
Definition at line 126 of file FP420SD.h.
Referenced by SetToLocal().
G4ThreeVector FP420SD::theExitPoint [private] |
Definition at line 127 of file FP420SD.h.
Referenced by SetToLocalExit().
FP420G4HitCollection* FP420SD::theHC [private] |
Definition at line 139 of file FP420SD.h.
Referenced by EndOfEvent(), HitExists(), Initialize(), PrintAll(), ProcessHits(), and StoreHit().
const SimTrackManager* FP420SD::theManager [private] |
float FP420SD::ThetaAtEntry [private] |
Definition at line 167 of file FP420SD.h.
Referenced by CreateNewHit(), and GetStepInfo().
G4Track* FP420SD::theTrack [private] |
Definition at line 144 of file FP420SD.h.
Referenced by CreateNewHit(), and GetStepInfo().
float FP420SD::Tof [private] |
Definition at line 163 of file FP420SD.h.
Referenced by CreateNewHit(), and GetStepInfo().
G4int FP420SD::tsID [private] |
Definition at line 142 of file FP420SD.h.
Referenced by CreateNewHit(), HitExists(), Initialize(), and UpdateHit().
G4double FP420SD::tSlice [private] |
Definition at line 151 of file FP420SD.h.
Referenced by CreateNewHit(), and GetStepInfo().
G4int FP420SD::tSliceID [private] |
Definition at line 148 of file FP420SD.h.
Referenced by CreateNewHit(), GetStepInfo(), HitExists(), and UpdateHit().
uint32_t FP420SD::unitID [private] |
Definition at line 147 of file FP420SD.h.
Referenced by CreateNewHit(), GetStepInfo(), HitExists(), and UpdateHit().
float FP420SD::Vx [private] |
Definition at line 171 of file FP420SD.h.
Referenced by CreateNewHit(), and GetStepInfo().
float FP420SD::Vy [private] |
Definition at line 171 of file FP420SD.h.
Referenced by CreateNewHit(), and GetStepInfo().
float FP420SD::Vz [private] |
Definition at line 171 of file FP420SD.h.
Referenced by CreateNewHit(), and GetStepInfo().
float FP420SD::X [private] |
Definition at line 172 of file FP420SD.h.
Referenced by CreateNewHit(), and GetStepInfo().
float FP420SD::Y [private] |
Definition at line 172 of file FP420SD.h.
Referenced by CreateNewHit(), and GetStepInfo().
float FP420SD::Z [private] |
Definition at line 172 of file FP420SD.h.
Referenced by CreateNewHit(), and GetStepInfo().