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