CMS 3D CMS Logo

CMSEmStandardPhysicsXS.cc
Go to the documentation of this file.
3 #include "G4EmParameters.hh"
4 #include "G4ParticleTable.hh"
5 
6 #include "G4ParticleDefinition.hh"
7 #include "G4LossTableManager.hh"
8 #include "G4EmParameters.hh"
9 
10 #include "G4ComptonScattering.hh"
11 #include "G4GammaConversion.hh"
12 #include "G4PhotoElectricEffect.hh"
13 #include "G4LivermorePhotoElectricModel.hh"
14 #include "G4KleinNishinaModel.hh"
15 
16 #include "G4hMultipleScattering.hh"
17 #include "G4eMultipleScattering.hh"
18 #include "G4MuMultipleScattering.hh"
19 #include "G4CoulombScattering.hh"
20 #include "G4eCoulombScatteringModel.hh"
21 #include "G4WentzelVIModel.hh"
22 #include "G4UrbanMscModel.hh"
23 #include "G4GoudsmitSaundersonMscModel.hh"
24 #include "G4MscStepLimitType.hh"
25 
26 #include "G4eIonisation.hh"
27 #include "G4eBremsstrahlung.hh"
28 #include "G4eplusAnnihilation.hh"
29 #include "G4Generator2BS.hh"
30 #include "G4SeltzerBergerModel.hh"
31 
32 #include "G4MuIonisation.hh"
33 #include "G4MuBremsstrahlung.hh"
34 #include "G4MuPairProduction.hh"
35 
36 #include "G4MuBremsstrahlungModel.hh"
37 #include "G4MuPairProductionModel.hh"
38 #include "G4hBremsstrahlungModel.hh"
39 #include "G4hPairProductionModel.hh"
40 #include "G4ePairProduction.hh"
41 
42 #include "G4hIonisation.hh"
43 #include "G4ionIonisation.hh"
44 #include "G4hBremsstrahlung.hh"
45 #include "G4hPairProduction.hh"
46 #include "G4UAtomicDeexcitation.hh"
47 
48 #include "G4NuclearStopping.hh"
49 
50 #include "G4Gamma.hh"
51 #include "G4Electron.hh"
52 #include "G4Positron.hh"
53 #include "G4MuonPlus.hh"
54 #include "G4MuonMinus.hh"
55 #include "G4TauMinus.hh"
56 #include "G4TauPlus.hh"
57 #include "G4PionPlus.hh"
58 #include "G4PionMinus.hh"
59 #include "G4KaonPlus.hh"
60 #include "G4KaonMinus.hh"
61 #include "G4BMesonMinus.hh"
62 #include "G4BMesonPlus.hh"
63 #include "G4DMesonMinus.hh"
64 #include "G4DMesonPlus.hh"
65 #include "G4Proton.hh"
66 #include "G4AntiProton.hh"
67 #include "G4SigmaMinus.hh"
68 #include "G4AntiSigmaMinus.hh"
69 #include "G4SigmaPlus.hh"
70 #include "G4AntiSigmaPlus.hh"
71 #include "G4XiMinus.hh"
72 #include "G4AntiXiMinus.hh"
73 #include "G4OmegaMinus.hh"
74 #include "G4AntiOmegaMinus.hh"
75 #include "G4LambdacPlus.hh"
76 #include "G4AntiLambdacPlus.hh"
77 #include "G4XicPlus.hh"
78 #include "G4AntiXicPlus.hh"
79 #include "G4Deuteron.hh"
80 #include "G4Triton.hh"
81 #include "G4He3.hh"
82 #include "G4Alpha.hh"
83 #include "G4GenericIon.hh"
84 
85 #include "G4PhysicsListHelper.hh"
86 #include "G4BuilderType.hh"
87 #include "G4RegionStore.hh"
88 #include "G4Region.hh"
89 
90 #include "G4SystemOfUnits.hh"
91 
92 CMSEmStandardPhysicsXS::CMSEmStandardPhysicsXS(G4int ver) : G4VPhysicsConstructor("CMSEmStandard_emn"), verbose(ver) {
93  G4EmParameters* param = G4EmParameters::Instance();
94  param->SetDefaults();
95  param->SetVerbose(verbose);
96  param->SetApplyCuts(true);
97  param->SetLowestElectronEnergy(100 * eV);
98  param->SetStepFunction(0.8, 1 * CLHEP::mm);
99  param->SetUseMottCorrection(true); // use Mott-correction for e-/e+ msc gs
100  param->SetMscRangeFactor(0.2);
101  param->SetMscStepLimitType(fMinimal);
102  param->SetFluo(true);
103  SetPhysicsType(bElectromagnetic);
104 }
105 
107 
109  // gamma
110  G4Gamma::Gamma();
111 
112  // leptons
114  G4Positron::Positron();
115  G4MuonPlus::MuonPlus();
116  G4MuonMinus::MuonMinus();
117  G4TauMinus::TauMinusDefinition();
118  G4TauPlus::TauPlusDefinition();
119 
120  // mesons
121  G4PionPlus::PionPlusDefinition();
122  G4PionMinus::PionMinusDefinition();
123  G4KaonPlus::KaonPlusDefinition();
124  G4KaonMinus::KaonMinusDefinition();
125  G4DMesonMinus::DMesonMinusDefinition();
126  G4DMesonPlus::DMesonPlusDefinition();
127  G4BMesonMinus::BMesonMinusDefinition();
128  G4BMesonPlus::BMesonPlusDefinition();
129 
130  // barions
131  G4Proton::Proton();
132  G4AntiProton::AntiProton();
133  G4SigmaMinus::SigmaMinusDefinition();
134  G4AntiSigmaMinus::AntiSigmaMinusDefinition();
135  G4SigmaPlus::SigmaPlusDefinition();
136  G4AntiSigmaPlus::AntiSigmaPlusDefinition();
137  G4XiMinus::XiMinusDefinition();
138  G4AntiXiMinus::AntiXiMinusDefinition();
139  G4OmegaMinus::OmegaMinusDefinition();
140  G4AntiOmegaMinus::AntiOmegaMinusDefinition();
141  G4LambdacPlus::LambdacPlusDefinition();
142  G4AntiLambdacPlus::AntiLambdacPlusDefinition();
143  G4XicPlus::XicPlusDefinition();
144  G4AntiXicPlus::AntiXicPlusDefinition();
145 
146  // ions
147  G4Deuteron::Deuteron();
148  G4Triton::Triton();
149  G4He3::He3();
150  G4Alpha::Alpha();
151  G4GenericIon::GenericIonDefinition();
152 }
153 
155  if (verbose > 0) {
156  G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
157  }
158 
159  // This EM builder takes default models of Geant4 10 EMV.
160  // Multiple scattering by Urban for all particles
161  // except e+e- below 100 MeV for which the Urban model is used
162 
163  G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
164 
165  // muon & hadron bremsstrahlung and pair production
166  G4MuBremsstrahlung* mub = nullptr;
167  G4MuPairProduction* mup = nullptr;
168  G4hBremsstrahlung* pib = nullptr;
169  G4hPairProduction* pip = nullptr;
170  G4hBremsstrahlung* kb = nullptr;
171  G4hPairProduction* kp = nullptr;
172  G4hBremsstrahlung* pb = nullptr;
173  G4hPairProduction* pp = nullptr;
174  G4ePairProduction* ee = nullptr;
175 
176  // muon & hadron multiple scattering
177  G4MuMultipleScattering* mumsc = nullptr;
178  G4hMultipleScattering* pimsc = nullptr;
179  G4hMultipleScattering* kmsc = nullptr;
180  G4hMultipleScattering* hmsc = nullptr;
181 
182  // muon and hadron single scattering
183  G4CoulombScattering* muss = nullptr;
184  G4CoulombScattering* piss = nullptr;
185  G4CoulombScattering* kss = nullptr;
186 
187  // high energy limit for e+- scattering models and bremsstrahlung
188  G4double highEnergyLimit = 100 * MeV;
189 
190  // nuclear stopping
191  G4NuclearStopping* pnuc = nullptr;
192 
193  G4Region* aRegion = G4RegionStore::GetInstance()->GetRegion("HcalRegion", false);
194  G4Region* bRegion = G4RegionStore::GetInstance()->GetRegion("HGCalRegion", false);
195 
196  G4ParticleTable* table = G4ParticleTable::GetParticleTable();
197  EmParticleList emList;
198  for (const auto& particleName : emList.PartNames()) {
199  G4ParticleDefinition* particle = table->FindParticle(particleName);
200 
201  if (particleName == "gamma") {
202  G4PhotoElectricEffect* photo = new G4PhotoElectricEffect();
203  photo->SetEmModel(new G4LivermorePhotoElectricModel());
204  ph->RegisterProcess(photo, particle);
205  G4ComptonScattering* compt = new G4ComptonScattering();
206  compt->SetEmModel(new G4KleinNishinaModel());
207  ph->RegisterProcess(compt, particle);
208  ph->RegisterProcess(new G4GammaConversion(), particle);
209 
210  } else if (particleName == "e-") {
211  G4eIonisation* eioni = new G4eIonisation();
212 
213  G4eMultipleScattering* msc = new G4eMultipleScattering;
214  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
215  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
216  G4GoudsmitSaundersonMscModel* msc3 = new G4GoudsmitSaundersonMscModel();
217  msc3->SetStepLimitType(fUseSafetyPlus);
218  msc3->SetRangeFactor(0.08);
219  msc3->SetSkin(3.0);
220  msc1->SetHighEnergyLimit(highEnergyLimit);
221  msc2->SetLowEnergyLimit(highEnergyLimit);
222  msc3->SetHighEnergyLimit(highEnergyLimit);
223  msc3->SetLocked(true);
224  msc->SetEmModel(msc1);
225  msc->SetEmModel(msc2);
226  msc->AddEmModel(-1, msc3, aRegion);
227  if (bRegion)
228  msc->AddEmModel(-1, msc3, bRegion);
229 
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  // bremsstrahlung
238  G4eBremsstrahlung* brem = new G4eBremsstrahlung();
239  G4SeltzerBergerModel* br1 = new G4SeltzerBergerModel();
240  G4eBremsstrahlungRelModel* br2 = new G4eBremsstrahlungRelModel();
241  br1->SetAngularDistribution(new G4Generator2BS());
242  br2->SetAngularDistribution(new G4Generator2BS());
243  brem->SetEmModel(br1);
244  brem->SetEmModel(br2);
245  br1->SetHighEnergyLimit(GeV);
246 
247  if (!ee) {
248  ee = new G4ePairProduction();
249  }
250 
251  ph->RegisterProcess(msc, particle);
252  ph->RegisterProcess(eioni, particle);
253  ph->RegisterProcess(brem, particle);
254  ph->RegisterProcess(ee, particle);
255  ph->RegisterProcess(ss, particle);
256 
257  } else if (particleName == "e+") {
258  G4eIonisation* eioni = new G4eIonisation();
259 
260  G4eMultipleScattering* msc = new G4eMultipleScattering;
261  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
262  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
263  G4GoudsmitSaundersonMscModel* msc3 = new G4GoudsmitSaundersonMscModel();
264  msc3->SetStepLimitType(fUseSafetyPlus);
265  msc3->SetRangeFactor(0.08);
266  msc3->SetSkin(3.0);
267  msc1->SetHighEnergyLimit(highEnergyLimit);
268  msc2->SetLowEnergyLimit(highEnergyLimit);
269  msc3->SetHighEnergyLimit(highEnergyLimit);
270  msc3->SetLocked(true);
271  msc->SetEmModel(msc1);
272  msc->SetEmModel(msc2);
273  msc->AddEmModel(-1, msc3, aRegion);
274  if (bRegion)
275  msc->AddEmModel(-1, msc3, bRegion);
276 
277  G4eCoulombScatteringModel* ssm = new G4eCoulombScatteringModel();
278  G4CoulombScattering* ss = new G4CoulombScattering();
279  ss->SetEmModel(ssm);
280  ss->SetMinKinEnergy(highEnergyLimit);
281  ssm->SetLowEnergyLimit(highEnergyLimit);
282  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
283 
284  // bremsstrahlung
285  G4eBremsstrahlung* brem = new G4eBremsstrahlung();
286  G4SeltzerBergerModel* br1 = new G4SeltzerBergerModel();
287  G4eBremsstrahlungRelModel* br2 = new G4eBremsstrahlungRelModel();
288  br1->SetAngularDistribution(new G4Generator2BS());
289  br2->SetAngularDistribution(new G4Generator2BS());
290  brem->SetEmModel(br1);
291  brem->SetEmModel(br2);
292  br1->SetHighEnergyLimit(GeV);
293 
294  if (!ee) {
295  ee = new G4ePairProduction();
296  }
297 
298  ph->RegisterProcess(msc, particle);
299  ph->RegisterProcess(eioni, particle);
300  ph->RegisterProcess(brem, particle);
301  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
302  ph->RegisterProcess(ee, particle);
303  ph->RegisterProcess(ss, particle);
304 
305  } else if (particleName == "mu+" || particleName == "mu-") {
306  if (nullptr == mub) {
307  mub = new G4MuBremsstrahlung();
308  mup = new G4MuPairProduction();
309  mumsc = new G4MuMultipleScattering();
310  mumsc->SetEmModel(new G4WentzelVIModel());
311  muss = new G4CoulombScattering();
312  }
313  ph->RegisterProcess(mumsc, particle);
314  ph->RegisterProcess(new G4MuIonisation(), particle);
315  ph->RegisterProcess(mub, particle);
316  ph->RegisterProcess(mup, particle);
317  ph->RegisterProcess(muss, particle);
318 
319  } else if (particleName == "alpha" || particleName == "He3") {
320  if (!pnuc) {
321  pnuc = new G4NuclearStopping();
322  }
323 
324  ph->RegisterProcess(new G4hMultipleScattering(), particle);
325  ph->RegisterProcess(new G4ionIonisation(), particle);
326  ph->RegisterProcess(pnuc, particle);
327 
328  } else if (particleName == "GenericIon") {
329  if (nullptr == hmsc) {
330  hmsc = new G4hMultipleScattering("ionmsc");
331  }
332  if (!pnuc) {
333  pnuc = new G4NuclearStopping();
334  }
335  ph->RegisterProcess(hmsc, particle);
336  ph->RegisterProcess(new G4ionIonisation(), particle);
337  ph->RegisterProcess(pnuc, particle);
338 
339  } else if (particleName == "pi+" || particleName == "pi-") {
340  if (nullptr == pib) {
341  pib = new G4hBremsstrahlung();
342  pip = new G4hPairProduction();
343  pimsc = new G4hMultipleScattering();
344  pimsc->SetEmModel(new G4WentzelVIModel());
345  piss = new G4CoulombScattering();
346  }
347  ph->RegisterProcess(pimsc, particle);
348  ph->RegisterProcess(new G4hIonisation(), particle);
349  ph->RegisterProcess(pib, particle);
350  ph->RegisterProcess(pip, particle);
351  ph->RegisterProcess(piss, particle);
352 
353  } else if (particleName == "kaon+" || particleName == "kaon-") {
354  if (nullptr == kb) {
355  kb = new G4hBremsstrahlung();
356  kp = new G4hPairProduction();
357  kmsc = new G4hMultipleScattering();
358  kmsc->SetEmModel(new G4WentzelVIModel());
359  kss = new G4CoulombScattering();
360  }
361  ph->RegisterProcess(kmsc, particle);
362  ph->RegisterProcess(new G4hIonisation(), particle);
363  ph->RegisterProcess(kb, particle);
364  ph->RegisterProcess(kp, particle);
365  ph->RegisterProcess(kss, particle);
366 
367  } else if (particleName == "proton" || particleName == "anti_proton") {
368  if (nullptr == pb) {
369  pb = new G4hBremsstrahlung();
370  pp = new G4hPairProduction();
371  }
372  if (!pnuc) {
373  pnuc = new G4NuclearStopping();
374  }
375 
376  G4hMultipleScattering* pmsc = new G4hMultipleScattering();
377  pmsc->SetEmModel(new G4WentzelVIModel());
378  G4hIonisation* hIoni = new G4hIonisation();
379  G4CoulombScattering* pss = new G4CoulombScattering();
380 
381  ph->RegisterProcess(pmsc, particle);
382  ph->RegisterProcess(hIoni, particle);
383  ph->RegisterProcess(pb, particle);
384  ph->RegisterProcess(pp, particle);
385  ph->RegisterProcess(pss, particle);
386 
387  } else if (particleName == "B+" || particleName == "B-" || particleName == "D+" || particleName == "D-" ||
388  particleName == "Ds+" || particleName == "Ds-" || particleName == "anti_He3" ||
389  particleName == "anti_alpha" || particleName == "anti_deuteron" || particleName == "anti_lambda_c+" ||
390  particleName == "anti_omega-" || particleName == "anti_sigma_c+" || particleName == "anti_sigma_c++" ||
391  particleName == "anti_sigma+" || particleName == "anti_sigma-" || particleName == "anti_triton" ||
392  particleName == "anti_xi_c+" || particleName == "anti_xi-" || particleName == "deuteron" ||
393  particleName == "lambda_c+" || particleName == "omega-" || particleName == "sigma_c+" ||
394  particleName == "sigma_c++" || particleName == "sigma+" || particleName == "sigma-" ||
395  particleName == "tau+" || particleName == "tau-" || particleName == "triton" ||
396  particleName == "xi_c+" || particleName == "xi-") {
397  if (nullptr == hmsc) {
398  hmsc = new G4hMultipleScattering("ionmsc");
399  }
400  ph->RegisterProcess(hmsc, particle);
401  ph->RegisterProcess(new G4hIonisation(), particle);
402  }
403  }
404  if (pnuc) {
405  pnuc->SetMaxKinEnergy(MeV);
406  }
407  // Deexcitation
408  //
409  G4VAtomDeexcitation* de = new G4UAtomicDeexcitation();
410  G4LossTableManager::Instance()->SetAtomDeexcitation(de);
411 }
const double GeV
Definition: MathUtil.h:16
const std::vector< G4String > & PartNames() const
std::pair< ALIstring, ALIstring > pss
Definition: Fit.h:27
int kp
susybsm::MuonSegmentCollection msc
Definition: classes.h:32
const double MeV
dbl * Gamma
Definition: mlp_gen.cc:38