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