CMS 3D CMS Logo

CMSEmStandardPhysicsXS.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 "G4PhysicsListHelper.hh"
48 #include "G4BuilderType.hh"
49 #include "G4GammaGeneralProcess.hh"
50 
51 #include "G4Version.hh"
52 #if G4VERSION_NUMBER >= 1110
53 #include "G4ProcessManager.hh"
54 #include "G4TransportationWithMsc.hh"
55 #endif
56 
57 #include "G4RegionStore.hh"
58 #include "G4Region.hh"
59 #include "G4GammaGeneralProcess.hh"
60 
61 #include "G4SystemOfUnits.hh"
62 
64  : G4VPhysicsConstructor("CMSEmStandard_emn") {
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->SetApplyCuts(true);
71  param->SetMinEnergy(100 * CLHEP::eV);
72  param->SetNumberOfBinsPerDecade(20);
73  param->SetStepFunction(0.8, 1 * CLHEP::mm);
74  param->SetMscRangeFactor(0.2);
75  param->SetMscStepLimitType(fMinimal);
76  param->SetFluo(true);
77  param->SetUseMottCorrection(true); // use Mott-correction for e-/e+ msc gs
78  SetPhysicsType(bElectromagnetic);
79  fRangeFactor = p.getParameter<double>("G4MscRangeFactor");
80  fGeomFactor = p.getParameter<double>("G4MscGeomFactor");
81  fSafetyFactor = p.getParameter<double>("G4MscSafetyFactor");
82  fLambdaLimit = p.getParameter<double>("G4MscLambdaLimit") * CLHEP::mm;
83  std::string msc = p.getParameter<std::string>("G4MscStepLimit");
84  fStepLimitType = fUseSafety;
85  if (msc == "UseSafetyPlus") {
86  fStepLimitType = fUseSafetyPlus;
87  }
88  if (msc == "Minimal") {
89  fStepLimitType = fMinimal;
90  }
91  double tcut = p.getParameter<double>("G4TrackingCut") * CLHEP::MeV;
92  param->SetLowestElectronEnergy(tcut);
93  param->SetLowestMuHadEnergy(tcut);
94 }
95 
97  // minimal set of particles for EM physics
98  G4EmBuilder::ConstructMinimalEmSet();
99 }
100 
102  if (verboseLevel > 0) {
103  edm::LogVerbatim("PhysicsList") << "### " << GetPhysicsName() << " Construct Processes";
104  }
105 
106  // This EM builder takes default models of Geant4 10 EMV.
107  // Multiple scattering by Urban for all particles
108  // except e+e- below 100 MeV for which the Urban model is used
109  G4EmBuilder::PrepareEMPhysics();
110 
111  G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
112 
113  // processes used by several particles
114  G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
115  G4NuclearStopping* pnuc(nullptr);
116 
117  // high energy limit for e+- scattering models
118  G4double highEnergyLimit = G4EmParameters::Instance()->MscEnergyLimit();
119  const G4Region* aRegion = G4RegionStore::GetInstance()->GetRegion("HcalRegion", false);
120  const G4Region* bRegion = G4RegionStore::GetInstance()->GetRegion("HGCalRegion", false);
121 
122  // Add gamma EM Processes
123  G4ParticleDefinition* particle = G4Gamma::Gamma();
124 
125  // Photoelectric
126  G4PhotoElectricEffect* pe = new G4PhotoElectricEffect();
127  G4VEmModel* theLivermorePEModel = new G4LivermorePhotoElectricModel();
128  pe->SetEmModel(theLivermorePEModel);
129 
130  // Compton scattering
131  G4ComptonScattering* cs = new G4ComptonScattering;
132  cs->SetEmModel(new G4KleinNishinaModel());
133 
134  // Gamma conversion
135  G4GammaConversion* gc = new G4GammaConversion();
136  G4VEmModel* conv = new G4BetheHeitler5DModel();
137  gc->SetEmModel(conv);
138 
139  if (G4EmParameters::Instance()->GeneralProcessActive()) {
141  sp->AddEmProcess(pe);
142  sp->AddEmProcess(cs);
143  sp->AddEmProcess(gc);
144  sp->AddEmProcess(new G4RayleighScattering());
145  G4LossTableManager::Instance()->SetGammaGeneralProcess(sp);
146  ph->RegisterProcess(sp, particle);
147  } else {
148  ph->RegisterProcess(pe, particle);
149  ph->RegisterProcess(cs, particle);
150  ph->RegisterProcess(gc, particle);
151  ph->RegisterProcess(new G4RayleighScattering(), particle);
152  }
153 
154  // e-
155  particle = G4Electron::Electron();
156 
157  // multiple scattering
158  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
159  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
160  msc1->SetHighEnergyLimit(highEnergyLimit);
161  msc2->SetLowEnergyLimit(highEnergyLimit);
162 
163  // msc for HCAL using the Urban model
164  G4UrbanMscModel* msc4 = nullptr;
165  if (nullptr != aRegion) {
166  msc4 = new G4UrbanMscModel();
167  msc4->SetHighEnergyLimit(highEnergyLimit);
168  msc4->SetRangeFactor(fRangeFactor);
169  msc4->SetGeomFactor(fGeomFactor);
170  msc4->SetSafetyFactor(fSafetyFactor);
171  msc4->SetLambdaLimit(fLambdaLimit);
172  msc4->SetStepLimitType(fStepLimitType);
173  msc4->SetLocked(true);
174  }
175 
176  // msc GS with Mott-correction
177  G4GoudsmitSaundersonMscModel* msc3 = nullptr;
178  if (nullptr != bRegion) {
179  msc3 = new G4GoudsmitSaundersonMscModel();
180  msc3->SetHighEnergyLimit(highEnergyLimit);
181  msc3->SetRangeFactor(0.08);
182  msc3->SetSkin(3.);
183  msc3->SetStepLimitType(fUseSafetyPlus);
184  msc3->SetLocked(true);
185  }
186 
187 #if G4VERSION_NUMBER >= 1110
188  G4TransportationWithMscType transportationWithMsc = G4EmParameters::Instance()->TransportationWithMsc();
189  if (transportationWithMsc != G4TransportationWithMscType::fDisabled) {
190  G4ProcessManager* procManager = particle->GetProcessManager();
191  // Remove default G4Transportation and replace with G4TransportationWithMsc.
192  G4VProcess* removed = procManager->RemoveProcess(0);
193  if (removed->GetProcessName() != "Transportation") {
194  G4Exception("CMSEmStandardPhysics::ConstructProcess",
195  "em0050",
196  FatalException,
197  "replaced process is not G4Transportation!");
198  }
199  G4TransportationWithMsc* transportWithMsc =
201  if (transportationWithMsc == G4TransportationWithMscType::fMultipleSteps) {
202  transportWithMsc->SetMultipleSteps(true);
203  }
204  transportWithMsc->AddMscModel(msc1);
205  transportWithMsc->AddMscModel(msc2);
206  if (msc4 != nullptr) {
207  transportWithMsc->AddMscModel(msc4, -1, aRegion);
208  }
209  if (msc3 != nullptr) {
210  transportWithMsc->AddMscModel(msc3, -1, bRegion);
211  }
212  procManager->AddProcess(transportWithMsc, -1, 0, 0);
213  } else
214 #endif
215  {
216  // Register as a separate process.
217  G4eMultipleScattering* msc = new G4eMultipleScattering;
218  msc->SetEmModel(msc1);
219  msc->SetEmModel(msc2);
220  if (msc4 != nullptr) {
221  msc->AddEmModel(-1, msc4, aRegion);
222  }
223  if (msc3 != nullptr) {
224  msc->AddEmModel(-1, msc3, bRegion);
225  }
226  ph->RegisterProcess(msc, particle);
227  }
228 
229  // single scattering
230  G4eCoulombScatteringModel* ssm = new G4eCoulombScatteringModel();
231  G4CoulombScattering* ss = new G4CoulombScattering();
232  ss->SetEmModel(ssm);
233  ss->SetMinKinEnergy(highEnergyLimit);
234  ssm->SetLowEnergyLimit(highEnergyLimit);
235  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
236 
237  // ionisation
238  G4eIonisation* eioni = new G4eIonisation();
239 
240  // bremsstrahlung
241  G4eBremsstrahlung* brem = new G4eBremsstrahlung();
242  G4SeltzerBergerModel* br1 = new G4SeltzerBergerModel();
243  G4eBremsstrahlungRelModel* br2 = new G4eBremsstrahlungRelModel();
244  br1->SetAngularDistribution(new G4Generator2BS());
245  br2->SetAngularDistribution(new G4Generator2BS());
246  brem->SetEmModel(br1);
247  brem->SetEmModel(br2);
248  br1->SetHighEnergyLimit(CLHEP::GeV);
249 
250  G4ePairProduction* ee = new G4ePairProduction();
251 
252  // register processes
253  ph->RegisterProcess(eioni, particle);
254  ph->RegisterProcess(brem, particle);
255  ph->RegisterProcess(ee, particle);
256  ph->RegisterProcess(ss, particle);
257 
258  // e+
259  particle = G4Positron::Positron();
260 
261  // multiple scattering
262  msc1 = new G4UrbanMscModel();
263  msc2 = new G4WentzelVIModel();
264  msc1->SetHighEnergyLimit(highEnergyLimit);
265  msc2->SetLowEnergyLimit(highEnergyLimit);
266 
267  // msc for HCAL using the Urban model
268  if (nullptr != aRegion) {
269  msc4 = new G4UrbanMscModel();
270  msc4->SetHighEnergyLimit(highEnergyLimit);
271  msc4->SetRangeFactor(fRangeFactor);
272  msc4->SetGeomFactor(fGeomFactor);
273  msc4->SetSafetyFactor(fSafetyFactor);
274  msc4->SetLambdaLimit(fLambdaLimit);
275  msc4->SetStepLimitType(fStepLimitType);
276  msc4->SetLocked(true);
277  }
278 
279  // msc GS with Mott-correction
280  if (nullptr != bRegion) {
281  msc3 = new G4GoudsmitSaundersonMscModel();
282  msc3->SetHighEnergyLimit(highEnergyLimit);
283  msc3->SetRangeFactor(0.08);
284  msc3->SetSkin(3.);
285  msc3->SetStepLimitType(fUseSafetyPlus);
286  msc3->SetLocked(true);
287  }
288 
289 #if G4VERSION_NUMBER >= 1110
290  if (transportationWithMsc != G4TransportationWithMscType::fDisabled) {
291  G4ProcessManager* procManager = particle->GetProcessManager();
292  // Remove default G4Transportation and replace with G4TransportationWithMsc.
293  G4VProcess* removed = procManager->RemoveProcess(0);
294  if (removed->GetProcessName() != "Transportation") {
295  G4Exception("CMSEmStandardPhysics::ConstructProcess",
296  "em0050",
297  FatalException,
298  "replaced process is not G4Transportation!");
299  }
300  G4TransportationWithMsc* transportWithMsc =
302  if (transportationWithMsc == G4TransportationWithMscType::fMultipleSteps) {
303  transportWithMsc->SetMultipleSteps(true);
304  }
305  transportWithMsc->AddMscModel(msc1);
306  transportWithMsc->AddMscModel(msc2);
307  if (msc4 != nullptr) {
308  transportWithMsc->AddMscModel(msc4, -1, aRegion);
309  }
310  if (msc3 != nullptr) {
311  transportWithMsc->AddMscModel(msc3, -1, bRegion);
312  }
313  procManager->AddProcess(transportWithMsc, -1, 0, 0);
314  } else
315 #endif
316  {
317  // Register as a separate process.
318  G4eMultipleScattering* msc = new G4eMultipleScattering;
319  msc->SetEmModel(msc1);
320  msc->SetEmModel(msc2);
321  if (msc4 != nullptr) {
322  msc->AddEmModel(-1, msc4, aRegion);
323  }
324  if (msc3 != nullptr) {
325  msc->AddEmModel(-1, msc3, bRegion);
326  }
327  ph->RegisterProcess(msc, particle);
328  }
329 
330  // single scattering
331  ssm = new G4eCoulombScatteringModel();
332  ss = new G4CoulombScattering();
333  ss->SetEmModel(ssm);
334  ss->SetMinKinEnergy(highEnergyLimit);
335  ssm->SetLowEnergyLimit(highEnergyLimit);
336  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
337 
338  // ionisation
339  eioni = new G4eIonisation();
340 
341  // bremsstrahlung
342  brem = new G4eBremsstrahlung();
343  br1 = new G4SeltzerBergerModel();
344  br2 = new G4eBremsstrahlungRelModel();
345  br1->SetAngularDistribution(new G4Generator2BS());
346  br2->SetAngularDistribution(new G4Generator2BS());
347  brem->SetEmModel(br1);
348  brem->SetEmModel(br2);
349  br1->SetHighEnergyLimit(CLHEP::GeV);
350 
351  // register processes
352  ph->RegisterProcess(eioni, particle);
353  ph->RegisterProcess(brem, particle);
354  ph->RegisterProcess(ee, particle);
355  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
356  ph->RegisterProcess(ss, particle);
357 
358  // generic ion
359  particle = G4GenericIon::GenericIon();
360  G4ionIonisation* ionIoni = new G4ionIonisation();
361  ph->RegisterProcess(hmsc, particle);
362  ph->RegisterProcess(ionIoni, particle);
363 
364  // muons, hadrons, ions
365  G4EmBuilder::ConstructCharged(hmsc, pnuc);
366 }
Log< level::Info, true > LogVerbatim
CMSEmStandardPhysicsXS(G4int ver, const edm::ParameterSet &p)
EPOS::IO_EPOS conv
G4MscStepLimitType fStepLimitType