CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ParametrisedEMPhysics.cc
Go to the documentation of this file.
1 //
2 // Joanna Weng 08.2005
3 // Physics process for Gflash parameterisation
4 // modified by Soon Yung Jun, Dongwook Jang
5 // V.Ivanchenko rename the class, cleanup, and move
6 // to SimG4Core/Application - 2012/08/14
7 
12 
13 #include "G4FastSimulationManagerProcess.hh"
14 #include "G4ProcessManager.hh"
15 
16 #include "G4LeptonConstructor.hh"
17 #include "G4MesonConstructor.hh"
18 #include "G4BaryonConstructor.hh"
19 #include "G4ShortLivedConstructor.hh"
20 #include "G4IonConstructor.hh"
21 #include "G4RegionStore.hh"
22 #include "G4Electron.hh"
23 #include "G4Positron.hh"
24 
25 #include "G4EmProcessOptions.hh"
26 #include "G4PhysicsListHelper.hh"
27 #include "G4SystemOfUnits.hh"
28 #include "G4UAtomicDeexcitation.hh"
29 #include "G4LossTableManager.hh"
30 
32  : G4VPhysicsConstructor(name), theParSet(p)
33 {
34  theEMShowerModel = 0;
36 }
37 
39  delete theEMShowerModel;
40  delete theHadShowerModel;
41 }
42 
44 {
45  G4LeptonConstructor pLeptonConstructor;
46  pLeptonConstructor.ConstructParticle();
47 
48  G4MesonConstructor pMesonConstructor;
49  pMesonConstructor.ConstructParticle();
50 
51  G4BaryonConstructor pBaryonConstructor;
52  pBaryonConstructor.ConstructParticle();
53 
54  G4ShortLivedConstructor pShortLivedConstructor;
55  pShortLivedConstructor.ConstructParticle();
56 
57  G4IonConstructor pConstructor;
58  pConstructor.ConstructParticle();
59 }
60 
62 
63  // GFlash part
64  bool gem = theParSet.getParameter<bool>("GflashEcal");
65  bool ghad = theParSet.getParameter<bool>("GflashHcal");
66 
67  if(gem || ghad) {
68  edm::LogInfo("SimG4CoreApplication")
69  << "ParametrisedEMPhysics: GFlash Construct: "
70  << gem << " " << ghad;
71  G4FastSimulationManagerProcess * theFastSimulationManagerProcess =
72  new G4FastSimulationManagerProcess();
73  aParticleIterator->reset();
74  while ((*aParticleIterator)()) {
75  G4ParticleDefinition * particle = aParticleIterator->value();
76  G4ProcessManager * pmanager = particle->GetProcessManager();
77  G4String pname = particle->GetParticleName();
78  if(pname == "e-" || pname == "e+") {
79  pmanager->AddProcess(theFastSimulationManagerProcess, -1, -1, 1);
80  }
81  }
82 
83  if(gem) {
84  G4Region* aRegion =
85  G4RegionStore::GetInstance()->GetRegion("EcalRegion");
86 
87  if(!aRegion){
88  edm::LogInfo("SimG4CoreApplication")
89  << "ParametrisedEMPhysics::ConstructProcess: "
90  << "EcalRegion is not defined, GFlash will not be enabled for ECAL!";
91 
92  } else {
93 
94  //Electromagnetic Shower Model for ECAL
96  new GFlashEMShowerModel("GflashEMShowerModel",aRegion,theParSet);
97  //std::cout << "GFlash is defined for EcalRegion" << std::endl;
98  }
99  }
100  if(ghad) {
101  G4Region* aRegion =
102  G4RegionStore::GetInstance()->GetRegion("HcalRegion");
103  if(!aRegion) {
104  edm::LogInfo("SimG4CoreApplication")
105  << "ParametrisedEMPhysics::ConstructProcess: "
106  << "HcalRegion is not defined, GFlash will not be enabled for HCAL!";
107 
108  } else {
109 
110  //Electromagnetic Shower Model for HCAL
112  new GFlashEMShowerModel("GflashHadShowerModel",aRegion,theParSet);
113  //std::cout << "GFlash is defined for HcalRegion" << std::endl;
114  }
115  }
116  }
117  // Russian Roulette part
118  G4EmProcessOptions opt;
119  const G4int NREG = 6;
120  const G4String rname[NREG] = {"EcalRegion", "HcalRegion", "MuonIron",
121  "PreshowerRegion","CastorRegion",
122  "DefaultRegionForTheWorld"};
123  G4double rrfact[NREG] = { 1.0 };
124 
125  // Russian roulette for e-
126  double energyLim =
127  theParSet.getParameter<double>("RusRoElectronEnergyLimit")*MeV;
128  if(energyLim > 0.0) {
129  rrfact[0] = theParSet.getParameter<double>("RusRoEcalElectron");
130  rrfact[1] = theParSet.getParameter<double>("RusRoHcalElectron");
131  rrfact[2] = theParSet.getParameter<double>("RusRoMuonIronElectron");
132  rrfact[3] = theParSet.getParameter<double>("RusRoPreShowerElectron");
133  rrfact[4] = theParSet.getParameter<double>("RusRoCastorElectron");
134  rrfact[5] = theParSet.getParameter<double>("RusRoWorldElectron");
135  for(int i=0; i<NREG; ++i) {
136  if(rrfact[i] < 1.0) {
137  opt.ActivateSecondaryBiasing("eIoni",rname[i],rrfact[i],energyLim);
138  opt.ActivateSecondaryBiasing("hIoni",rname[i],rrfact[i],energyLim);
139  //opt.ActivateSecondaryBiasing("muIoni",rname[i],rrfact[i],energyLim);
140  //opt.ActivateSecondaryBiasing("ionIoni",rname[i],rrfact[i],energyLim);
141  //opt.ActivateSecondaryBiasingForGamma("phot",rname[i],rrfact[i],energyLim);
142  //opt.ActivateSecondaryBiasingForGamma("compt",rname[i],rrfact[i],energyLim);
143  //opt.ActivateSecondaryBiasingForGamma("conv",rname[i],rrfact[i],energyLim);
144  edm::LogInfo("SimG4CoreApplication")
145  << "ParametrisedEMPhysics: Russian Roulette"
146  << " for e- Prob= " << rrfact[i]
147  << " Elimit(MeV)= " << energyLim/CLHEP::MeV
148  << " inside " << rname[i];
149  }
150  }
151  }
152 
153  // Step limiters for e+-
154  bool eLimiter = theParSet.getParameter<bool>("ElectronStepLimit");
155  bool rLimiter = theParSet.getParameter<bool>("ElectronRangeTest");
156  bool pLimiter = theParSet.getParameter<bool>("PositronStepLimit");
157 
158  if(eLimiter || rLimiter || pLimiter) {
159  G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
160 
161  if(eLimiter || rLimiter) {
165  ph->RegisterProcess(theElectronLimiter, G4Electron::Electron());
166  }
167 
168  if(pLimiter){
171  ph->RegisterProcess(theElectronLimiter, G4Positron::Positron());
172  }
173  }
174  // enable fluorescence
175  bool fluo = theParSet.getParameter<bool>("FlagFluo");
176  if(fluo) {
177  G4VAtomDeexcitation* de = new G4UAtomicDeexcitation();
178  G4LossTableManager::Instance()->SetAtomDeexcitation(de);
179  de->SetFluo(true);
180  }
181 }
ElectronLimiter * theElectronLimiter
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
void SetRangeCheckFlag(G4bool)
edm::ParameterSet theParSet
ElectronLimiter * thePositronLimiter
const double MeV
void SetFieldCheckFlag(G4bool)
GFlashEMShowerModel * theHadShowerModel
ParametrisedEMPhysics(std::string name, const edm::ParameterSet &p)
GFlashEMShowerModel * theEMShowerModel