CMS 3D CMS Logo

/data/doxygen/doxygen-1.7.3/gen/CMSSW_4_2_8/src/SimG4Core/GFlash/src/GflashEMShowerModel.cc

Go to the documentation of this file.
00001 //
00002 // initial setup : E.Barberio & Joanna Weng 
00003 // big changes : Soon Jun & Dongwook Jang
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   // Mininum energy cutoff to parameterize
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   // This will be changed accordingly when the way dealing with CaloRegion changes later.
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   // The parameterization starts inside crystals
00071   std::size_t pos1 = lv->GetName().find("EBRY");
00072   std::size_t pos2 = lv->GetName().find("EFRY");
00073   /*
00074   std::size_t pos3 = lv->GetName().find("HVQ");
00075   std::size_t pos4 = lv->GetName().find("HF");
00076   if(pos1 == std::string::npos && pos2 == std::string::npos &&
00077      pos3 == std::string::npos && pos4 == std::string::npos) return false;
00078   */
00079   //@@@for now, HF is not a part of Gflash Envelopes
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   // Kill the parameterised particle:
00091   fastStep.KillPrimaryTrack();
00092   fastStep.ProposePrimaryTrackPathLength(0.0);
00093 
00094   //input variables for GflashEMShowerProfile with showerType = 1,5 (shower starts inside crystals)
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   // Do actual parameterization. The result of parameterization is gflashHitList
00103   theProfile->initialize(showerType,energy,globalTime,charge,position,momentum);
00104   theProfile->parameterization();
00105 
00106   //make hits
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       //put touchable for each hit so that touchable history keeps track of each step.
00126       theGflashNavigator->LocateGlobalPointAndUpdateTouchableHandle(spotIter->getPosition(),G4ThreeVector(0,0,0),
00127                                                                     theGflashTouchableHandle, false);
00128       updateGflashStep(spotIter->getPosition(),spotIter->getTime());
00129 
00130       // if there is a watcher defined in a job and the flag is turned on
00131       if(theWatcherOn) {
00132         SteppingAction* userSteppingAction = (SteppingAction*) G4EventManager::GetEventManager()->GetUserSteppingAction();
00133         userSteppingAction->m_g4StepSignal(theGflashStep);
00134       }
00135 
00136       // Send G4Step information to Hit/Digi if the volume is sensitive
00137       // Copied from G4SteppingManager.cc
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   //exclude regions where geometry are complicated
00170   //+- one supermodule around the EB/EE boundary: 1.479 +- 0.0174*5 
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 G4int GflashEMShowerModel::findShowerType(const G4FastTrack& fastTrack)
00179 {
00180   // Initialization of longitudinal and lateral parameters for
00181   // hadronic showers. Simulation of the intrinsic fluctuations
00182 
00183   // type of hadron showers subject to the shower starting point (ssp)
00184   // showerType = -1 : default (invalid)
00185   // showerType =  0 : ssp before EBRY (barrel crystal)
00186   // showerType =  1 : ssp inside EBRY
00187   // showerType =  2 : ssp after  EBRY before HB
00188   // showerType =  3 : ssp inside HB
00189   // showerType =  4 : ssp before EFRY (endcap crystal)
00190   // showerType =  5 : ssp inside EFRY
00191   // showerType =  6 : ssp after  EFRY before HE
00192   // showerType =  7 : ssp inside HE
00193 
00194   G4TouchableHistory* touch = (G4TouchableHistory*)(fastTrack.GetPrimaryTrack()->GetTouchable());
00195   G4LogicalVolume* lv = touch->GetVolume()->GetLogicalVolume();
00196 
00197   std::size_t pos1  = lv->GetName().find("EBRY");
00198   std::size_t pos11 = lv->GetName().find("EWAL");
00199   std::size_t pos12 = lv->GetName().find("EWRA");
00200   std::size_t pos2  = lv->GetName().find("EFRY");
00201 
00202   G4ThreeVector position = fastTrack.GetPrimaryTrack()->GetPosition()/cm;
00203   Gflash::CalorimeterNumber kCalor = Gflash::getCalorimeterNumber(position);
00204 
00205   G4int showerType = -1;
00206 
00207   //central
00208   if (kCalor == Gflash::kESPM || kCalor == Gflash::kHB ) {
00209 
00210     G4double posRho = position.getRho();
00211 
00212     if(pos1 != std::string::npos || pos11 != std::string::npos || pos12 != std::string::npos ) {
00213       showerType = 1;
00214     }
00215     else {
00216       if(kCalor == Gflash::kESPM) {
00217         showerType = 2;
00218         if( posRho < Gflash::Rmin[Gflash::kESPM]+ Gflash::ROffCrystalEB ) showerType = 0;
00219       }
00220       else showerType = 3;
00221     }
00222 
00223   }
00224   //forward
00225   else if (kCalor == Gflash::kENCA || kCalor == Gflash::kHE) {
00226     if(pos2 != std::string::npos) {
00227       showerType = 5;
00228     }
00229     else {
00230       if(kCalor == Gflash::kENCA) {
00231         showerType = 6;
00232         if(fabs(position.getZ()) < Gflash::Zmin[Gflash::kENCA] + Gflash::ZOffCrystalEE) showerType = 4;
00233       }
00234       else showerType = 7;
00235     }
00236     //@@@need z-dependent correction on the mean energy reponse
00237   }
00238 
00239   return showerType;
00240 }
00241 */