CMS 3D CMS Logo

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