test
CMS 3D CMS Logo

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