CMS 3D CMS Logo

CMSEmStandardPhysicsLPM.cc
Go to the documentation of this file.
3 
5 
6 #include "G4EmParameters.hh"
7 #include "G4ParticleTable.hh"
8 
9 #include "G4ParticleDefinition.hh"
10 #include "G4LossTableManager.hh"
11 
12 #include "G4ComptonScattering.hh"
13 #include "G4GammaConversion.hh"
14 #include "G4PhotoElectricEffect.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 
24 #include "G4eIonisation.hh"
25 #include "G4eBremsstrahlung.hh"
26 #include "G4eplusAnnihilation.hh"
27 #include "G4UAtomicDeexcitation.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 "G4PionPlus.hh"
44 #include "G4PionMinus.hh"
45 #include "G4KaonPlus.hh"
46 #include "G4KaonMinus.hh"
47 #include "G4Proton.hh"
48 #include "G4AntiProton.hh"
49 #include "G4GenericIon.hh"
50 
51 #include "G4PhysicsListHelper.hh"
52 #include "G4BuilderType.hh"
53 #include "G4RegionStore.hh"
54 #include "G4Region.hh"
55 #include "G4GammaGeneralProcess.hh"
56 #include "G4EmBuilder.hh"
57 
58 #include "G4SystemOfUnits.hh"
59 
61  : G4VPhysicsConstructor("CMSEmStandard_emm") {
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->SetStepFunction(0.8, 1 * CLHEP::mm);
69  param->SetMscRangeFactor(0.2);
70  param->SetMscStepLimitType(fMinimal);
71  SetPhysicsType(bElectromagnetic);
72  double tcut = p.getParameter<double>("G4TrackingCut") * CLHEP::MeV;
73  param->SetLowestElectronEnergy(tcut);
74  param->SetLowestMuHadEnergy(tcut);
75 }
76 
78  // minimal set of particles for EM physics
79  G4EmBuilder::ConstructMinimalEmSet();
80 }
81 
83  if (verboseLevel > 1) {
84  edm::LogVerbatim("PhysicsList") << "### " << GetPhysicsName() << " Construct Processes";
85  }
86 
87  // This EM builder takes default models of Geant4 10 EMV.
88  // Multiple scattering by Urban for all particles
89  // except e+e- below 100 MeV for which the Urban model is used
90 
91  G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
92  G4LossTableManager* man = G4LossTableManager::Instance();
93 
94  // muon & hadron bremsstrahlung and pair production
95  G4MuBremsstrahlung* mub = nullptr;
96  G4MuPairProduction* mup = nullptr;
97  G4hBremsstrahlung* pib = nullptr;
98  G4hPairProduction* pip = nullptr;
99  G4hBremsstrahlung* kb = nullptr;
100  G4hPairProduction* kp = nullptr;
101  G4hBremsstrahlung* pb = nullptr;
102  G4hPairProduction* pp = nullptr;
103 
104  // muon & hadron multiple scattering
105  G4MuMultipleScattering* mumsc = nullptr;
106  G4hMultipleScattering* pimsc = nullptr;
107  G4hMultipleScattering* kmsc = nullptr;
108  G4hMultipleScattering* pmsc = nullptr;
109  G4hMultipleScattering* hmsc = nullptr;
110 
111  // muon and hadron single scattering
112  G4CoulombScattering* muss = nullptr;
113  G4CoulombScattering* piss = nullptr;
114  G4CoulombScattering* kss = nullptr;
115 
116  // high energy limit for e+- scattering models and bremsstrahlung
117  G4double highEnergyLimit = 100 * MeV;
118 
119  G4Region* aRegion = G4RegionStore::GetInstance()->GetRegion("HcalRegion", false);
120  G4Region* bRegion = G4RegionStore::GetInstance()->GetRegion("HGCalRegion", false);
121  if (verboseLevel > 1) {
122  edm::LogVerbatim("PhysicsList") << "CMSEmStandardPhysicsLPM: HcalRegion " << aRegion << "; HGCalRegion " << bRegion;
123  }
124  G4ParticleTable* table = G4ParticleTable::GetParticleTable();
125  EmParticleList emList;
126  for (const auto& particleName : emList.PartNames()) {
127  G4ParticleDefinition* particle = table->FindParticle(particleName);
128 
129  if (particleName == "gamma") {
130  G4PhotoElectricEffect* pee = new G4PhotoElectricEffect();
131 
132  if (G4EmParameters::Instance()->GeneralProcessActive()) {
134  sp->AddEmProcess(pee);
135  sp->AddEmProcess(new G4ComptonScattering());
136  sp->AddEmProcess(new G4GammaConversion());
137  man->SetGammaGeneralProcess(sp);
138  ph->RegisterProcess(sp, particle);
139 
140  } else {
141  ph->RegisterProcess(pee, particle);
142  ph->RegisterProcess(new G4ComptonScattering(), particle);
143  ph->RegisterProcess(new G4GammaConversion(), particle);
144  }
145 
146  } else if (particleName == "e-") {
147  G4eIonisation* eioni = new G4eIonisation();
148 
149  G4eMultipleScattering* msc = new G4eMultipleScattering;
150  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
151  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
152  msc1->SetHighEnergyLimit(highEnergyLimit);
153  msc2->SetLowEnergyLimit(highEnergyLimit);
154  msc->SetEmModel(msc1);
155  msc->SetEmModel(msc2);
156 
157  // e-/e+ msc for HCAL and HGCAL using the Urban model
158  if (nullptr != aRegion || nullptr != bRegion) {
159  G4UrbanMscModel* msc3 = new G4UrbanMscModel();
160  msc3->SetHighEnergyLimit(highEnergyLimit);
161  msc3->SetLocked(true);
162 
163  if (nullptr != aRegion) {
164  msc->AddEmModel(-1, msc3, aRegion);
165  }
166  if (nullptr != bRegion) {
167  msc->AddEmModel(-1, msc3, bRegion);
168  }
169  }
170 
171  G4eCoulombScatteringModel* ssm = new G4eCoulombScatteringModel();
172  G4CoulombScattering* ss = new G4CoulombScattering();
173  ss->SetEmModel(ssm);
174  ss->SetMinKinEnergy(highEnergyLimit);
175  ssm->SetLowEnergyLimit(highEnergyLimit);
176  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
177 
178  ph->RegisterProcess(msc, particle);
179  ph->RegisterProcess(eioni, particle);
180  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
181  ph->RegisterProcess(ss, particle);
182 
183  } else if (particleName == "e+") {
184  G4eIonisation* eioni = new G4eIonisation();
185 
186  G4eMultipleScattering* msc = new G4eMultipleScattering;
187  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
188  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
189  msc1->SetHighEnergyLimit(highEnergyLimit);
190  msc2->SetLowEnergyLimit(highEnergyLimit);
191  msc->SetEmModel(msc1);
192  msc->SetEmModel(msc2);
193 
194  // e-/e+ msc for HCAL and HGCAL using the Urban model
195  if (nullptr != aRegion || nullptr != bRegion) {
196  G4UrbanMscModel* msc3 = new G4UrbanMscModel();
197  msc3->SetHighEnergyLimit(highEnergyLimit);
198  msc3->SetLocked(true);
199 
200  if (nullptr != aRegion) {
201  msc->AddEmModel(-1, msc3, aRegion);
202  }
203  if (nullptr != bRegion) {
204  msc->AddEmModel(-1, msc3, bRegion);
205  }
206  }
207 
208  G4eCoulombScatteringModel* ssm = new G4eCoulombScatteringModel();
209  G4CoulombScattering* ss = new G4CoulombScattering();
210  ss->SetEmModel(ssm);
211  ss->SetMinKinEnergy(highEnergyLimit);
212  ssm->SetLowEnergyLimit(highEnergyLimit);
213  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
214 
215  ph->RegisterProcess(msc, particle);
216  ph->RegisterProcess(eioni, particle);
217  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
218  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
219  ph->RegisterProcess(ss, particle);
220 
221  } else if (particleName == "mu+" || particleName == "mu-") {
222  if (nullptr == mub) {
223  mub = new G4MuBremsstrahlung();
224  mup = new G4MuPairProduction();
225  mumsc = new G4MuMultipleScattering();
226  mumsc->SetEmModel(new G4WentzelVIModel());
227  muss = new G4CoulombScattering();
228  }
229  ph->RegisterProcess(mumsc, particle);
230  ph->RegisterProcess(new G4MuIonisation(), particle);
231  ph->RegisterProcess(mub, particle);
232  ph->RegisterProcess(mup, particle);
233  ph->RegisterProcess(muss, particle);
234 
235  } else if (particleName == "alpha" || particleName == "He3") {
236  ph->RegisterProcess(new G4hMultipleScattering(), particle);
237  ph->RegisterProcess(new G4ionIonisation(), particle);
238 
239  } else if (particleName == "GenericIon") {
240  if (nullptr == hmsc) {
241  hmsc = new G4hMultipleScattering("ionmsc");
242  }
243  ph->RegisterProcess(hmsc, particle);
244  ph->RegisterProcess(new G4ionIonisation(), particle);
245 
246  } else if (particleName == "pi+" || particleName == "pi-") {
247  if (nullptr == pib) {
248  pib = new G4hBremsstrahlung();
249  pip = new G4hPairProduction();
250  pimsc = new G4hMultipleScattering();
251  pimsc->SetEmModel(new G4WentzelVIModel());
252  piss = new G4CoulombScattering();
253  }
254  ph->RegisterProcess(pimsc, particle);
255  ph->RegisterProcess(new G4hIonisation(), particle);
256  ph->RegisterProcess(pib, particle);
257  ph->RegisterProcess(pip, particle);
258  ph->RegisterProcess(piss, particle);
259 
260  } else if (particleName == "kaon+" || particleName == "kaon-") {
261  if (nullptr == kb) {
262  kb = new G4hBremsstrahlung();
263  kp = new G4hPairProduction();
264  kmsc = new G4hMultipleScattering();
265  kmsc->SetEmModel(new G4WentzelVIModel());
266  kss = new G4CoulombScattering();
267  }
268  ph->RegisterProcess(kmsc, particle);
269  ph->RegisterProcess(new G4hIonisation(), particle);
270  ph->RegisterProcess(kb, particle);
271  ph->RegisterProcess(kp, particle);
272  ph->RegisterProcess(kss, particle);
273 
274  } else if (particleName == "proton" || particleName == "anti_proton") {
275  if (nullptr == pb) {
276  pb = new G4hBremsstrahlung();
277  pp = new G4hPairProduction();
278  }
279  pmsc = new G4hMultipleScattering();
280  pmsc->SetEmModel(new G4WentzelVIModel());
281 
282  ph->RegisterProcess(pmsc, particle);
283  ph->RegisterProcess(new G4hIonisation(), particle);
284  ph->RegisterProcess(pb, particle);
285  ph->RegisterProcess(pp, particle);
286  ph->RegisterProcess(new G4CoulombScattering(), particle);
287 
288  } else if (particle->GetPDGCharge() != 0.0) {
289  if (nullptr == hmsc) {
290  hmsc = new G4hMultipleScattering("ionmsc");
291  }
292  ph->RegisterProcess(hmsc, particle);
293  ph->RegisterProcess(new G4hIonisation(), particle);
294  }
295  }
296  edm::LogVerbatim("PhysicsList") << "CMSEmStandardPhysicsLPM: EM physics is instantiated";
297 }
Log< level::Info, true > LogVerbatim
CMSEmStandardPhysicsLPM(G4int ver, const edm::ParameterSet &p)
const std::vector< G4String > & PartNames() const