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