CMS 3D CMS Logo

CMSEmStandardPhysicsEMZ.cc
Go to the documentation of this file.
3 
4 #include "G4SystemOfUnits.hh"
5 #include "G4ParticleDefinition.hh"
6 #include "G4LossTableManager.hh"
7 #include "G4EmParameters.hh"
8 #include "G4EmBuilder.hh"
9 
10 #include "G4ComptonScattering.hh"
11 #include "G4GammaConversion.hh"
12 #include "G4PhotoElectricEffect.hh"
13 #include "G4RayleighScattering.hh"
14 #include "G4PEEffectFluoModel.hh"
15 #include "G4KleinNishinaModel.hh"
16 #include "G4LowEPComptonModel.hh"
17 #include "G4BetheHeitler5DModel.hh"
18 #include "G4LivermorePhotoElectricModel.hh"
19 
20 #include "G4eMultipleScattering.hh"
21 #include "G4hMultipleScattering.hh"
22 #include "G4MscStepLimitType.hh"
23 #include "G4UrbanMscModel.hh"
24 #include "G4GoudsmitSaundersonMscModel.hh"
25 #include "G4DummyModel.hh"
26 #include "G4WentzelVIModel.hh"
27 #include "G4CoulombScattering.hh"
28 #include "G4eCoulombScatteringModel.hh"
29 
30 #include "G4eIonisation.hh"
31 #include "G4eBremsstrahlung.hh"
32 #include "G4Generator2BS.hh"
33 #include "G4SeltzerBergerModel.hh"
34 #include "G4ePairProduction.hh"
35 #include "G4UniversalFluctuation.hh"
36 
37 #include "G4eplusAnnihilation.hh"
38 
39 #include "G4hIonisation.hh"
40 #include "G4ionIonisation.hh"
41 
42 #include "G4Gamma.hh"
43 #include "G4Electron.hh"
44 #include "G4Positron.hh"
45 #include "G4GenericIon.hh"
46 
47 #include "G4IonParametrisedLossModel.hh"
48 #include "G4NuclearStopping.hh"
49 #include "G4ePairProduction.hh"
50 #include "G4LivermoreIonisationModel.hh"
51 #include "G4PenelopeIonisationModel.hh"
52 
53 #include "G4PhysicsListHelper.hh"
54 #include "G4BuilderType.hh"
55 #include "G4GammaGeneralProcess.hh"
56 
57 #include "G4RegionStore.hh"
58 #include "G4Region.hh"
59 #include "G4GammaGeneralProcess.hh"
60 
61 #include "G4SystemOfUnits.hh"
62 
63 CMSEmStandardPhysicsEMZ::CMSEmStandardPhysicsEMZ(G4int ver) : G4VPhysicsConstructor("CMSEmStandard_emz") {
64  SetVerboseLevel(ver);
65  G4EmParameters* param = G4EmParameters::Instance();
66  param->SetDefaults();
67  param->SetVerbose(ver);
68  param->SetMinEnergy(100 * CLHEP::eV);
69  param->SetLowestElectronEnergy(100 * CLHEP::eV);
70  param->SetNumberOfBinsPerDecade(20);
71  param->ActivateAngularGeneratorForIonisation(true);
72  param->SetStepFunction(0.2, 10 * CLHEP::um);
73  param->SetStepFunctionMuHad(0.1, 50 * CLHEP::um);
74  param->SetStepFunctionLightIons(0.1, 20 * CLHEP::um);
75  param->SetStepFunctionIons(0.1, 1 * CLHEP::um);
76  param->SetUseMottCorrection(true); // use Mott-correction for e-/e+ msc gs
77  param->SetMscStepLimitType(fUseSafetyPlus); // for e-/e+ msc gs
78  param->SetMscSkin(3); // error-free stepping for e-/e+ msc gs
79  param->SetMscRangeFactor(0.08); // error-free stepping for e-/e+ msc gs
80  param->SetMuHadLateralDisplacement(true);
81  param->SetFluo(true);
82  param->SetUseICRU90Data(true);
83  param->SetMaxNIELEnergy(1 * CLHEP::MeV);
84  SetPhysicsType(bElectromagnetic);
85 }
86 
88 
90  // minimal set of particles for EM physics
91  G4EmBuilder::ConstructMinimalEmSet();
92 }
93 
95  if (verboseLevel > 0) {
96  edm::LogVerbatim("PhysicsList") << "### " << GetPhysicsName() << " Construct Processes";
97  }
98 
99  // This EM builder takes default models of Geant4 10 EMV.
100  // Multiple scattering by Urban for all particles
101  // except e+e- below 100 MeV for which the Urban model is used
102  G4EmBuilder::PrepareEMPhysics();
103 
104  G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
105 
106  // processes used by several particles
107  G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
108  G4double nielEnergyLimit = G4EmParameters::Instance()->MaxNIELEnergy();
109  G4NuclearStopping* pnuc = nullptr;
110  if (nielEnergyLimit > 0.0) {
111  pnuc = new G4NuclearStopping();
112  pnuc->SetMaxKinEnergy(nielEnergyLimit);
113  }
114 
115  // high energy limit for e+- scattering models
116  G4double highEnergyLimit = G4EmParameters::Instance()->MscEnergyLimit();
117 
118  // Add gamma EM Processes
119  G4ParticleDefinition* particle = G4Gamma::Gamma();
120 
121  // Photoelectric
122  G4PhotoElectricEffect* pe = new G4PhotoElectricEffect();
123  G4VEmModel* theLivermorePEModel = new G4LivermorePhotoElectricModel();
124  pe->SetEmModel(theLivermorePEModel);
125 
126  // Compton scattering
127  G4ComptonScattering* cs = new G4ComptonScattering;
128  cs->SetEmModel(new G4KleinNishinaModel());
129  G4VEmModel* theLowEPComptonModel = new G4LowEPComptonModel();
130  theLowEPComptonModel->SetHighEnergyLimit(20 * CLHEP::MeV);
131  cs->AddEmModel(0, theLowEPComptonModel);
132 
133  // Gamma conversion
134  G4GammaConversion* gc = new G4GammaConversion();
135  G4VEmModel* conv = new G4BetheHeitler5DModel();
136  gc->SetEmModel(conv);
137 
138  if (G4EmParameters::Instance()->GeneralProcessActive()) {
139  G4GammaGeneralProcess* sp = new G4GammaGeneralProcess();
140  sp->AddEmProcess(pe);
141  sp->AddEmProcess(cs);
142  sp->AddEmProcess(gc);
143  sp->AddEmProcess(new G4RayleighScattering());
144  G4LossTableManager::Instance()->SetGammaGeneralProcess(sp);
145  ph->RegisterProcess(sp, particle);
146  } else {
147  ph->RegisterProcess(pe, particle);
148  ph->RegisterProcess(cs, particle);
149  ph->RegisterProcess(gc, particle);
150  ph->RegisterProcess(new G4RayleighScattering(), particle);
151  }
152 
153  // e-
154  particle = G4Electron::Electron();
155 
156  // multiple scattering
157  G4eMultipleScattering* msc = new G4eMultipleScattering();
158  // e-/e+ msc gs with Mott-correction
159  // (Mott-correction is set through G4EmParameters)
160  G4GoudsmitSaundersonMscModel* msc1 = new G4GoudsmitSaundersonMscModel();
161  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
162  msc1->SetHighEnergyLimit(highEnergyLimit);
163  msc2->SetLowEnergyLimit(highEnergyLimit);
164  msc->SetEmModel(msc1);
165  msc->SetEmModel(msc2);
166 
167  G4eCoulombScatteringModel* ssm = new G4eCoulombScatteringModel();
168  G4CoulombScattering* ss = new G4CoulombScattering();
169  ss->SetEmModel(ssm);
170  ss->SetMinKinEnergy(highEnergyLimit);
171  ssm->SetLowEnergyLimit(highEnergyLimit);
172  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
173 
174  // single scattering
175  ssm = new G4eCoulombScatteringModel();
176  ss = new G4CoulombScattering();
177  ss->SetEmModel(ssm);
178  ss->SetMinKinEnergy(highEnergyLimit);
179  ssm->SetLowEnergyLimit(highEnergyLimit);
180  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
181 
182  // ionisation
183  G4eIonisation* eioni = new G4eIonisation();
184  G4VEmModel* theIoniLiv = new G4LivermoreIonisationModel();
185  theIoniLiv->SetHighEnergyLimit(0.1 * CLHEP::MeV);
186  eioni->AddEmModel(0, theIoniLiv, new G4UniversalFluctuation());
187 
188  // bremsstrahlung
189  G4eBremsstrahlung* brem = new G4eBremsstrahlung();
190  G4SeltzerBergerModel* br1 = new G4SeltzerBergerModel();
191  G4eBremsstrahlungRelModel* br2 = new G4eBremsstrahlungRelModel();
192  br1->SetAngularDistribution(new G4Generator2BS());
193  br2->SetAngularDistribution(new G4Generator2BS());
194  brem->SetEmModel(br1);
195  brem->SetEmModel(br2);
196  br1->SetHighEnergyLimit(CLHEP::GeV);
197 
198  G4ePairProduction* ee = new G4ePairProduction();
199 
200  // register processes
201  ph->RegisterProcess(msc, particle);
202  ph->RegisterProcess(eioni, particle);
203  ph->RegisterProcess(brem, particle);
204  ph->RegisterProcess(ee, particle);
205  ph->RegisterProcess(ss, particle);
206 
207  // e+
208  particle = G4Positron::Positron();
209 
210  // multiple scattering
211  msc = new G4eMultipleScattering();
212  // e-/e+ msc gs with Mott-correction
213  // (Mott-correction is set through G4EmParameters)
214  msc1 = new G4GoudsmitSaundersonMscModel();
215  msc2 = new G4WentzelVIModel();
216  msc1->SetHighEnergyLimit(highEnergyLimit);
217  msc2->SetLowEnergyLimit(highEnergyLimit);
218  msc->SetEmModel(msc1);
219  msc->SetEmModel(msc2);
220 
221  // single scattering
222  ssm = new G4eCoulombScatteringModel();
223  ss = new G4CoulombScattering();
224  ss->SetEmModel(ssm);
225  ss->SetMinKinEnergy(highEnergyLimit);
226  ssm->SetLowEnergyLimit(highEnergyLimit);
227  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
228 
229  // ionisation
230  eioni = new G4eIonisation();
231  /*
232  G4VEmModel* pen = new G4PenelopeIonisationModel();
233  pen->SetHighEnergyLimit(0.1*CLHEP::MeV);
234  eioni->AddEmModel(0, pen, new G4UniversalFluctuation());
235  */
236  // bremsstrahlung
237  brem = new G4eBremsstrahlung();
238  br1 = new G4SeltzerBergerModel();
239  br2 = new G4eBremsstrahlungRelModel();
240  br1->SetAngularDistribution(new G4Generator2BS());
241  br2->SetAngularDistribution(new G4Generator2BS());
242  brem->SetEmModel(br1);
243  brem->SetEmModel(br2);
244  br1->SetHighEnergyLimit(CLHEP::GeV);
245 
246  // register processes
247  ph->RegisterProcess(msc, particle);
248  ph->RegisterProcess(eioni, particle);
249  ph->RegisterProcess(brem, particle);
250  ph->RegisterProcess(ee, particle);
251  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
252  ph->RegisterProcess(ss, particle);
253 
254  // generic ion
255  particle = G4GenericIon::GenericIon();
256  G4ionIonisation* ionIoni = new G4ionIonisation();
257  ionIoni->SetEmModel(new G4IonParametrisedLossModel());
258  ph->RegisterProcess(hmsc, particle);
259  ph->RegisterProcess(ionIoni, particle);
260  if (nullptr != pnuc)
261  ph->RegisterProcess(pnuc, particle);
262 
263  // muons, hadrons, ions
264  G4EmBuilder::ConstructCharged(hmsc, pnuc);
265 }
Log< level::Info, true > LogVerbatim
EPOS::IO_EPOS conv