![]() |
![]() |
#include <GFlashEMShowerModel.h>
Public Member Functions | |
void | DoIt (const G4FastTrack &, G4FastStep &) |
GFlashEMShowerModel (const G4String &name, G4Envelope *env, edm::ParameterSet parSet) | |
G4bool | IsApplicable (const G4ParticleDefinition &) |
G4bool | ModelTrigger (const G4FastTrack &) |
virtual | ~GFlashEMShowerModel () |
Private Member Functions | |
G4bool | excludeDetectorRegion (const G4FastTrack &fastTrack) |
void | makeHits (const G4FastTrack &fastTrack) |
void | updateGflashStep (G4ThreeVector position, G4double time) |
Private Attributes | |
G4Navigator * | theGflashNavigator |
G4Step * | theGflashStep |
G4TouchableHandle | theGflashTouchableHandle |
edm::ParameterSet | theParSet |
GflashEMShowerProfile * | theProfile |
const G4Region * | theRegion |
bool | theWatcherOn |
Definition at line 29 of file GFlashEMShowerModel.h.
GFlashEMShowerModel::GFlashEMShowerModel | ( | const G4String & | name, |
G4Envelope * | env, | ||
edm::ParameterSet | parSet | ||
) |
Definition at line 24 of file GFlashEMShowerModel.cc.
References edm::ParameterSet::getParameter(), theGflashNavigator, theGflashStep, theGflashTouchableHandle, theProfile, theRegion, and theWatcherOn.
: G4VFastSimulationModel(modelName, envelope), theParSet(parSet) { theWatcherOn = parSet.getParameter<bool>("watcherOn"); theProfile = new GflashEMShowerProfile(parSet); theRegion = const_cast<const G4Region*>(envelope); theGflashStep = new G4Step(); theGflashTouchableHandle = new G4TouchableHistory(); theGflashNavigator = new G4Navigator(); }
GFlashEMShowerModel::~GFlashEMShowerModel | ( | ) | [virtual] |
Definition at line 42 of file GFlashEMShowerModel.cc.
References theGflashStep, and theProfile.
{ delete theProfile; delete theGflashStep; }
void GFlashEMShowerModel::DoIt | ( | const G4FastTrack & | fastTrack, |
G4FastStep & | fastStep | ||
) |
Definition at line 88 of file GFlashEMShowerModel.cc.
References DeDxDiscriminatorTools::charge(), relval_parameters_module::energy, Gflash::findShowerType(), GflashEMShowerProfile::initialize(), makeHits(), GflashEMShowerProfile::parameterization(), position, and theProfile.
{ // Kill the parameterised particle: fastStep.KillPrimaryTrack(); fastStep.ProposePrimaryTrackPathLength(0.0); // Input variables for GFlashEMShowerProfile with showerType = 1,5 // Shower starts inside crystals G4double energy = fastTrack.GetPrimaryTrack()->GetKineticEnergy()/GeV; G4double globalTime = fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetGlobalTime(); G4double charge = fastTrack.GetPrimaryTrack()->GetStep()->GetPreStepPoint()->GetCharge(); G4ThreeVector position = fastTrack.GetPrimaryTrack()->GetPosition()/cm; G4ThreeVector momentum = fastTrack.GetPrimaryTrack()->GetMomentum()/GeV; G4int showerType = Gflash::findShowerType(position); // Do actual parameterization // The result of parameterization is gflashHitList theProfile->initialize(showerType,energy,globalTime,charge, position,momentum); theProfile->parameterization(); // Make hits makeHits(fastTrack); }
G4bool GFlashEMShowerModel::excludeDetectorRegion | ( | const G4FastTrack & | fastTrack | ) | [private] |
Definition at line 175 of file GFlashEMShowerModel.cc.
References eta.
Referenced by ModelTrigger().
{ G4bool isExcluded=false; //exclude regions where geometry are complicated //+- one supermodule around the EB/EE boundary: 1.479 +- 0.0174*5 G4double eta = fastTrack.GetPrimaryTrack()->GetPosition().pseudoRapidity(); if(std::fabs(eta) > 1.392 && std::fabs(eta) < 1.566) { return true; } return isExcluded; }
G4bool GFlashEMShowerModel::IsApplicable | ( | const G4ParticleDefinition & | particleType | ) |
Definition at line 49 of file GFlashEMShowerModel.cc.
{
return ( &particleType == G4Electron::Electron() ||
&particleType == G4Positron::Positron() );
}
void GFlashEMShowerModel::makeHits | ( | const G4FastTrack & | fastTrack | ) | [private] |
Definition at line 117 of file GFlashEMShowerModel.cc.
References GflashEMShowerProfile::getGflashHitList(), SteppingAction::m_g4StepSignal, theGflashNavigator, theGflashStep, theGflashTouchableHandle, theProfile, theRegion, theWatcherOn, and updateGflashStep().
Referenced by DoIt().
{ std::vector<GflashHit>& gflashHitList = theProfile->getGflashHitList(); theGflashStep->SetTrack(const_cast<G4Track*>(fastTrack.GetPrimaryTrack())); theGflashStep->GetPostStepPoint() ->SetProcessDefinedStep(const_cast<G4VProcess*>(fastTrack.GetPrimaryTrack() ->GetStep()->GetPostStepPoint()->GetProcessDefinedStep())); theGflashNavigator->SetWorldVolume(G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume()); std::vector<GflashHit>::const_iterator spotIter = gflashHitList.begin(); std::vector<GflashHit>::const_iterator spotIterEnd = gflashHitList.end(); for(; spotIter != spotIterEnd; spotIter++){ // Put touchable for each hit so that touchable history // keeps track of each step. theGflashNavigator->LocateGlobalPointAndUpdateTouchableHandle(spotIter->getPosition(),G4ThreeVector(0,0,0),theGflashTouchableHandle, false); updateGflashStep(spotIter->getPosition(),spotIter->getTime()); // If there is a watcher defined in a job and the flag is turned on if(theWatcherOn) { SteppingAction* userSteppingAction = (SteppingAction*) G4EventManager::GetEventManager()->GetUserSteppingAction(); userSteppingAction->m_g4StepSignal(theGflashStep); } // Send G4Step information to Hit/Digi if the volume is sensitive // Copied from G4SteppingManager.cc G4VPhysicalVolume* aCurrentVolume = theGflashStep->GetPreStepPoint()->GetPhysicalVolume(); if( aCurrentVolume == 0 ) { continue; } G4LogicalVolume* lv = aCurrentVolume->GetLogicalVolume(); if(lv->GetRegion() != theRegion) { continue; } theGflashStep->GetPreStepPoint()->SetSensitiveDetector(aCurrentVolume->GetLogicalVolume()->GetSensitiveDetector()); G4VSensitiveDetector* aSensitive = theGflashStep->GetPreStepPoint()->GetSensitiveDetector(); if( aSensitive == 0 ) { continue; } theGflashStep->SetTotalEnergyDeposit(spotIter->getEnergy()); aSensitive->Hit(theGflashStep); } }
G4bool GFlashEMShowerModel::ModelTrigger | ( | const G4FastTrack & | fastTrack | ) |
Definition at line 56 of file GFlashEMShowerModel.cc.
References excludeDetectorRegion(), and theRegion.
{ // Mininum energy cutoff to parameterize if(fastTrack.GetPrimaryTrack()->GetKineticEnergy() < GeV) { return false; } if(excludeDetectorRegion(fastTrack)) { return false; } // This will be changed accordingly when the way // dealing with CaloRegion changes later. G4TouchableHistory* touch = (G4TouchableHistory*)(fastTrack.GetPrimaryTrack()->GetTouchable()); G4VPhysicalVolume* pCurrentVolume = touch->GetVolume(); if( pCurrentVolume == 0) { return false; } G4LogicalVolume* lv = pCurrentVolume->GetLogicalVolume(); if(lv->GetRegion() != theRegion) { return false; } //std::cout << "GFlashEMShowerModel::ModelTrigger: LV " // << lv->GetRegion()->GetName() << std::endl; // The parameterization starts inside crystals //std::size_t pos1 = lv->GetName().find("EBRY"); //std::size_t pos2 = lv->GetName().find("EFRY"); //std::size_t pos3 = lv->GetName().find("HVQ"); //std::size_t pos4 = lv->GetName().find("HF"); //if(pos1 == std::string::npos && pos2 == std::string::npos && // pos3 == std::string::npos && pos4 == std::string::npos) return false; return true; }
void GFlashEMShowerModel::updateGflashStep | ( | G4ThreeVector | position, |
G4double | time | ||
) | [private] |
Definition at line 165 of file GFlashEMShowerModel.cc.
References theGflashStep, and theGflashTouchableHandle.
Referenced by makeHits().
{ theGflashStep->GetPostStepPoint()->SetGlobalTime(timeGlobal); theGflashStep->GetPreStepPoint()->SetPosition(spotPosition); theGflashStep->GetPostStepPoint()->SetPosition(spotPosition); theGflashStep->GetPreStepPoint()->SetTouchableHandle(theGflashTouchableHandle); }
G4Navigator* GFlashEMShowerModel::theGflashNavigator [private] |
Definition at line 56 of file GFlashEMShowerModel.h.
Referenced by GFlashEMShowerModel(), and makeHits().
G4Step* GFlashEMShowerModel::theGflashStep [private] |
Definition at line 55 of file GFlashEMShowerModel.h.
Referenced by GFlashEMShowerModel(), makeHits(), updateGflashStep(), and ~GFlashEMShowerModel().
G4TouchableHandle GFlashEMShowerModel::theGflashTouchableHandle [private] |
Definition at line 57 of file GFlashEMShowerModel.h.
Referenced by GFlashEMShowerModel(), makeHits(), and updateGflashStep().
Definition at line 48 of file GFlashEMShowerModel.h.
Definition at line 51 of file GFlashEMShowerModel.h.
Referenced by DoIt(), GFlashEMShowerModel(), makeHits(), and ~GFlashEMShowerModel().
const G4Region* GFlashEMShowerModel::theRegion [private] |
Definition at line 53 of file GFlashEMShowerModel.h.
Referenced by GFlashEMShowerModel(), makeHits(), and ModelTrigger().
bool GFlashEMShowerModel::theWatcherOn [private] |
Definition at line 49 of file GFlashEMShowerModel.h.
Referenced by GFlashEMShowerModel(), and makeHits().