11 #include "G4Electron.hh"
12 #include "G4Positron.hh"
13 #include "G4VProcess.hh"
14 #include "G4VPhysicalVolume.hh"
15 #include "G4LogicalVolume.hh"
16 #include "G4TransportationManager.hh"
17 #include "G4EventManager.hh"
18 #include "G4FastSimulationManager.hh"
19 #include "G4TouchableHandle.hh"
20 #include "G4VSensitiveDetector.hh"
23 : G4VFastSimulationModel(modelName, envelope), theParSet(parSet) {
45 return ( &particleType == G4Electron::ElectronDefinition() ||
46 &particleType == G4Positron::PositronDefinition() );
54 if(fastTrack.GetPrimaryTrack()->GetKineticEnergy() < 1.0*GeV)
return false;
57 G4bool trigger = fastTrack.GetPrimaryTrack()->GetDefinition() == G4Electron::ElectronDefinition() ||
58 fastTrack.GetPrimaryTrack()->GetDefinition() == G4Positron::PositronDefinition();
60 if(!trigger)
return false;
63 G4TouchableHistory* touch = (G4TouchableHistory*)(fastTrack.GetPrimaryTrack()->GetTouchable());
64 G4VPhysicalVolume* pCurrentVolume = touch->GetVolume();
65 if( pCurrentVolume == 0)
return false;
67 G4LogicalVolume* lv = pCurrentVolume->GetLogicalVolume();
68 if(lv->GetRegion()->GetName() !=
"CaloRegion")
return false;
71 std::size_t pos1 = lv->GetName().find(
"EBRY");
72 std::size_t pos2 = lv->GetName().find(
"EFRY");
80 if(pos1 == std::string::npos && pos2 == std::string::npos )
return false;
91 fastStep.KillPrimaryTrack();
92 fastStep.ProposePrimaryTrackPathLength(0.0);
95 G4double
energy = fastTrack.GetPrimaryTrack()->GetKineticEnergy()/GeV;
96 G4double globalTime = fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetGlobalTime();
97 G4double
charge = fastTrack.GetPrimaryTrack()->GetStep()->GetPreStepPoint()->GetCharge();
98 G4ThreeVector
position = fastTrack.GetPrimaryTrack()->GetPosition() / cm;
99 G4ThreeVector momentum = fastTrack.GetPrimaryTrack()->GetMomentum()/GeV;
114 theGflashStep->SetTrack(const_cast<G4Track*>(fastTrack.GetPrimaryTrack()));
116 theGflashStep->GetPostStepPoint()->SetProcessDefinedStep(const_cast<G4VProcess*>
117 (fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetProcessDefinedStep()));
118 theGflashNavigator->SetWorldVolume(G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume());
120 std::vector<GflashHit>::const_iterator spotIter = gflashHitList.begin();
121 std::vector<GflashHit>::const_iterator spotIterEnd = gflashHitList.end();
123 for( ; spotIter != spotIterEnd; spotIter++){
126 theGflashNavigator->LocateGlobalPointAndUpdateTouchableHandle(spotIter->getPosition(),G4ThreeVector(0,0,0),
139 G4VPhysicalVolume* aCurrentVolume =
theGflashStep->GetPreStepPoint()->GetPhysicalVolume();
140 if( aCurrentVolume == 0 )
continue;
142 G4LogicalVolume* lv = aCurrentVolume->GetLogicalVolume();
143 if(lv->GetRegion()->GetName() !=
"CaloRegion")
continue;
145 theGflashStep->GetPreStepPoint()->SetSensitiveDetector(aCurrentVolume->GetLogicalVolume()->GetSensitiveDetector());
146 G4VSensitiveDetector* aSensitive =
theGflashStep->GetPreStepPoint()->GetSensitiveDetector();
148 if( aSensitive == 0 )
continue;
158 theGflashStep->GetPostStepPoint()->SetGlobalTime(timeGlobal);
160 theGflashStep->GetPostStepPoint()->SetPosition(spotPosition);
167 G4bool isExcluded=
false;
171 G4double
eta = fastTrack.GetPrimaryTrack()->GetPosition().pseudoRapidity() ;
172 if(std::fabs(eta) > 1.392 && std::fabs(eta) < 1.566)
return true;
T getParameter(std::string const &) const
std::vector< GflashHit > & getGflashHitList()
SimActivityRegistry::G4StepSignal m_g4StepSignal
void DoIt(const G4FastTrack &, G4FastStep &)
G4bool excludeDetectorRegion(const G4FastTrack &fastTrack)
void updateGflashStep(G4ThreeVector position, G4double time)
G4TouchableHandle theGflashTouchableHandle
static int position[TOTALCHAMBERS][3]
int findShowerType(const Gflash3Vector position)
GflashEMShowerProfile * theProfile
void initialize(int showerType, double energy, double globalTime, double charge, Gflash3Vector &position, Gflash3Vector &momentum)
G4bool IsApplicable(const G4ParticleDefinition &)
G4Navigator * theGflashNavigator
void makeHits(const G4FastTrack &fastTrack)
GflashEMShowerModel(G4String name, G4Envelope *env, edm::ParameterSet parSet)
G4bool ModelTrigger(const G4FastTrack &)