Go to the documentation of this file.00001
00002
00003
00004
00005 #include "SimG4Core/Application/interface/SteppingAction.h"
00006 #include "SimG4Core/GFlash/interface/GflashEMShowerModel.h"
00007
00008 #include "SimGeneral/GFlash/interface/GflashEMShowerProfile.h"
00009 #include "SimGeneral/GFlash/interface/GflashHit.h"
00010
00011 #include "G4Electron.hh"
00012 #include "G4Positron.hh"
00013 #include "G4VProcess.hh"
00014 #include "G4VPhysicalVolume.hh"
00015 #include "G4LogicalVolume.hh"
00016 #include "G4TransportationManager.hh"
00017 #include "G4EventManager.hh"
00018 #include "G4FastSimulationManager.hh"
00019 #include "G4TouchableHandle.hh"
00020 #include "G4VSensitiveDetector.hh"
00021
00022 GflashEMShowerModel::GflashEMShowerModel(G4String modelName, G4Envelope* envelope, edm::ParameterSet parSet)
00023 : G4VFastSimulationModel(modelName, envelope), theParSet(parSet) {
00024
00025 theWatcherOn = parSet.getParameter<bool>("watcherOn");
00026
00027 theProfile = new GflashEMShowerProfile(parSet);
00028
00029 theGflashStep = new G4Step();
00030 theGflashTouchableHandle = new G4TouchableHistory();
00031 theGflashNavigator = new G4Navigator();
00032
00033 }
00034
00035
00036
00037 GflashEMShowerModel::~GflashEMShowerModel() {
00038
00039 if(theProfile) delete theProfile;
00040 if(theGflashStep) delete theGflashStep;
00041 }
00042
00043 G4bool GflashEMShowerModel::IsApplicable(const G4ParticleDefinition& particleType) {
00044
00045 return ( &particleType == G4Electron::ElectronDefinition() ||
00046 &particleType == G4Positron::PositronDefinition() );
00047
00048 }
00049
00050
00051 G4bool GflashEMShowerModel::ModelTrigger(const G4FastTrack & fastTrack ) {
00052
00053
00054 if(fastTrack.GetPrimaryTrack()->GetKineticEnergy() < 1.0*GeV) return false;
00055 if(excludeDetectorRegion(fastTrack)) return false;
00056
00057 G4bool trigger = fastTrack.GetPrimaryTrack()->GetDefinition() == G4Electron::ElectronDefinition() ||
00058 fastTrack.GetPrimaryTrack()->GetDefinition() == G4Positron::PositronDefinition();
00059
00060 if(!trigger) return false;
00061
00062
00063 G4TouchableHistory* touch = (G4TouchableHistory*)(fastTrack.GetPrimaryTrack()->GetTouchable());
00064 G4VPhysicalVolume* pCurrentVolume = touch->GetVolume();
00065 if( pCurrentVolume == 0) return false;
00066
00067 G4LogicalVolume* lv = pCurrentVolume->GetLogicalVolume();
00068 if(lv->GetRegion()->GetName() != "CaloRegion") return false;
00069
00070
00071 std::size_t pos1 = lv->GetName().find("EBRY");
00072 std::size_t pos2 = lv->GetName().find("EFRY");
00073
00074
00075
00076
00077
00078
00079
00080 if(pos1 == std::string::npos && pos2 == std::string::npos ) return false;
00081
00082 return true;
00083
00084 }
00085
00086
00087
00088 void GflashEMShowerModel::DoIt(const G4FastTrack& fastTrack, G4FastStep& fastStep) {
00089
00090
00091 fastStep.KillPrimaryTrack();
00092 fastStep.ProposePrimaryTrackPathLength(0.0);
00093
00094
00095 G4double energy = fastTrack.GetPrimaryTrack()->GetKineticEnergy()/GeV;
00096 G4double globalTime = fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetGlobalTime();
00097 G4double charge = fastTrack.GetPrimaryTrack()->GetStep()->GetPreStepPoint()->GetCharge();
00098 G4ThreeVector position = fastTrack.GetPrimaryTrack()->GetPosition() / cm;
00099 G4ThreeVector momentum = fastTrack.GetPrimaryTrack()->GetMomentum()/GeV;
00100 G4int showerType = Gflash::findShowerType(position);
00101
00102
00103 theProfile->initialize(showerType,energy,globalTime,charge,position,momentum);
00104 theProfile->parameterization();
00105
00106
00107 makeHits(fastTrack);
00108 }
00109
00110 void GflashEMShowerModel::makeHits(const G4FastTrack& fastTrack) {
00111
00112 std::vector<GflashHit>& gflashHitList = theProfile->getGflashHitList();
00113
00114 theGflashStep->SetTrack(const_cast<G4Track*>(fastTrack.GetPrimaryTrack()));
00115
00116 theGflashStep->GetPostStepPoint()->SetProcessDefinedStep(const_cast<G4VProcess*>
00117 (fastTrack.GetPrimaryTrack()->GetStep()->GetPostStepPoint()->GetProcessDefinedStep()));
00118 theGflashNavigator->SetWorldVolume(G4TransportationManager::GetTransportationManager()->GetNavigatorForTracking()->GetWorldVolume());
00119
00120 std::vector<GflashHit>::const_iterator spotIter = gflashHitList.begin();
00121 std::vector<GflashHit>::const_iterator spotIterEnd = gflashHitList.end();
00122
00123 for( ; spotIter != spotIterEnd; spotIter++){
00124
00125
00126 theGflashNavigator->LocateGlobalPointAndUpdateTouchableHandle(spotIter->getPosition(),G4ThreeVector(0,0,0),
00127 theGflashTouchableHandle, false);
00128 updateGflashStep(spotIter->getPosition(),spotIter->getTime());
00129
00130
00131 if(theWatcherOn) {
00132 SteppingAction* userSteppingAction = (SteppingAction*) G4EventManager::GetEventManager()->GetUserSteppingAction();
00133 userSteppingAction->m_g4StepSignal(theGflashStep);
00134 }
00135
00136
00137
00138
00139 G4VPhysicalVolume* aCurrentVolume = theGflashStep->GetPreStepPoint()->GetPhysicalVolume();
00140 if( aCurrentVolume == 0 ) continue;
00141
00142 G4LogicalVolume* lv = aCurrentVolume->GetLogicalVolume();
00143 if(lv->GetRegion()->GetName() != "CaloRegion") continue;
00144
00145 theGflashStep->GetPreStepPoint()->SetSensitiveDetector(aCurrentVolume->GetLogicalVolume()->GetSensitiveDetector());
00146 G4VSensitiveDetector* aSensitive = theGflashStep->GetPreStepPoint()->GetSensitiveDetector();
00147
00148 if( aSensitive == 0 ) continue;
00149
00150 theGflashStep->SetTotalEnergyDeposit(spotIter->getEnergy());
00151 aSensitive->Hit(theGflashStep);
00152 }
00153
00154 }
00155
00156 void GflashEMShowerModel::updateGflashStep(G4ThreeVector spotPosition, G4double timeGlobal)
00157 {
00158 theGflashStep->GetPostStepPoint()->SetGlobalTime(timeGlobal);
00159 theGflashStep->GetPreStepPoint()->SetPosition(spotPosition);
00160 theGflashStep->GetPostStepPoint()->SetPosition(spotPosition);
00161 theGflashStep->GetPreStepPoint()->SetTouchableHandle(theGflashTouchableHandle);
00162 }
00163
00164
00165 G4bool GflashEMShowerModel::excludeDetectorRegion(const G4FastTrack& fastTrack) {
00166
00167 G4bool isExcluded=false;
00168
00169
00170
00171 G4double eta = fastTrack.GetPrimaryTrack()->GetPosition().pseudoRapidity() ;
00172 if(std::fabs(eta) > 1.392 && std::fabs(eta) < 1.566) return true;
00173
00174 return isExcluded;
00175 }
00176
00177
00178
00179
00180
00181
00182
00183
00184
00185
00186
00187
00188
00189
00190
00191
00192
00193
00194
00195
00196
00197
00198
00199
00200
00201
00202
00203
00204
00205
00206
00207
00208
00209
00210
00211
00212
00213
00214
00215
00216
00217
00218
00219
00220
00221
00222
00223
00224
00225
00226
00227
00228
00229
00230
00231
00232
00233
00234
00235
00236
00237
00238
00239
00240
00241