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