CMS 3D CMS Logo

CMSEmStandardPhysics.cc
Go to the documentation of this file.
3 
4 #include "G4SystemOfUnits.hh"
5 #include "G4ParticleDefinition.hh"
6 #include "G4EmParameters.hh"
7 #include "G4EmBuilder.hh"
8 
9 #include "G4ComptonScattering.hh"
10 #include "G4GammaConversion.hh"
11 #include "G4PhotoElectricEffect.hh"
12 #include "G4LivermorePhotoElectricModel.hh"
13 
14 #include "G4MscStepLimitType.hh"
15 
16 #include "G4eMultipleScattering.hh"
17 #include "G4hMultipleScattering.hh"
18 #include "G4eCoulombScatteringModel.hh"
19 #include "G4CoulombScattering.hh"
20 #include "G4WentzelVIModel.hh"
21 #include "G4UrbanMscModel.hh"
22 
23 #include "G4eIonisation.hh"
24 #include "G4eBremsstrahlung.hh"
25 #include "G4eplusAnnihilation.hh"
26 
27 #include "G4hIonisation.hh"
28 #include "G4ionIonisation.hh"
29 
30 #include "G4ParticleTable.hh"
31 #include "G4Gamma.hh"
32 #include "G4Electron.hh"
33 #include "G4Positron.hh"
34 #include "G4GenericIon.hh"
35 
36 #include "G4PhysicsListHelper.hh"
37 #include "G4BuilderType.hh"
38 #include "G4GammaGeneralProcess.hh"
39 #include "G4LossTableManager.hh"
40 
41 #include "G4Version.hh"
42 #if G4VERSION_NUMBER >= 1110
43 #include "G4ProcessManager.hh"
44 #include "G4TransportationWithMsc.hh"
45 #endif
46 
47 #include "G4RegionStore.hh"
48 #include "G4Region.hh"
49 #include <string>
50 
52  : G4VPhysicsConstructor("CMSEmStandard_emm") {
53  SetVerboseLevel(ver);
54  G4EmParameters* param = G4EmParameters::Instance();
55  param->SetDefaults();
56  param->SetVerbose(ver);
57  param->SetApplyCuts(true);
58  param->SetStepFunction(0.8, 1 * CLHEP::mm);
59  param->SetMscRangeFactor(0.2);
60  param->SetMscStepLimitType(fMinimal);
61  param->SetFluo(false);
62  SetPhysicsType(bElectromagnetic);
63  fRangeFactor = p.getParameter<double>("G4MscRangeFactor");
64  fGeomFactor = p.getParameter<double>("G4MscGeomFactor");
65  fSafetyFactor = p.getParameter<double>("G4MscSafetyFactor");
66  fLambdaLimit = p.getParameter<double>("G4MscLambdaLimit") * CLHEP::mm;
67  std::string msc = p.getParameter<std::string>("G4MscStepLimit");
68  fStepLimitType = fUseSafety;
69  if (msc == "UseSafetyPlus") {
70  fStepLimitType = fUseSafetyPlus;
71  }
72  if (msc == "Minimal") {
73  fStepLimitType = fMinimal;
74  }
75  double tcut = p.getParameter<double>("G4TrackingCut") * CLHEP::MeV;
76  param->SetLowestElectronEnergy(tcut);
77  param->SetLowestMuHadEnergy(tcut);
78 }
79 
81 
83  // minimal set of particles for EM physics
84  G4EmBuilder::ConstructMinimalEmSet();
85 }
86 
88  if (verboseLevel > 0) {
89  edm::LogVerbatim("PhysicsList") << "### " << GetPhysicsName() << " Construct EM Processes";
90  }
91 
92  // This EM builder takes default models of Geant4 10 EMV.
93  // Multiple scattering by WentzelVI for all particles except:
94  // a) e+e- below 100 MeV for which the Urban model is used
95  // b) ions for which Urban model is used
96  G4EmBuilder::PrepareEMPhysics();
97 
98  G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
99  // processes used by several particles
100  G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
101  G4NuclearStopping* pnuc(nullptr);
102 
103  // high energy limit for e+- scattering models
104  G4double highEnergyLimit = G4EmParameters::Instance()->MscEnergyLimit();
105 
106  const G4Region* aRegion = G4RegionStore::GetInstance()->GetRegion("HcalRegion", false);
107  const G4Region* bRegion = G4RegionStore::GetInstance()->GetRegion("HGCalRegion", false);
108 
109  // Add gamma EM Processes
110  G4ParticleDefinition* particle = G4Gamma::Gamma();
111 
112  G4PhotoElectricEffect* pee = new G4PhotoElectricEffect();
113 
114  if (G4EmParameters::Instance()->GeneralProcessActive()) {
116  sp->AddEmProcess(pee);
117  sp->AddEmProcess(new G4ComptonScattering());
118  sp->AddEmProcess(new G4GammaConversion());
119  G4LossTableManager::Instance()->SetGammaGeneralProcess(sp);
120  ph->RegisterProcess(sp, particle);
121 
122  } else {
123  ph->RegisterProcess(pee, particle);
124  ph->RegisterProcess(new G4ComptonScattering(), particle);
125  ph->RegisterProcess(new G4GammaConversion(), particle);
126  }
127 
128  // e-
129  particle = G4Electron::Electron();
130 
131  G4eIonisation* eioni = new G4eIonisation();
132 
133  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
134  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
135  msc1->SetHighEnergyLimit(highEnergyLimit);
136  msc2->SetLowEnergyLimit(highEnergyLimit);
137 
138  // e-/e+ msc for HCAL and HGCAL using the Urban model
139  G4UrbanMscModel* msc3 = nullptr;
140  if (nullptr != aRegion || nullptr != bRegion) {
141  msc3 = new G4UrbanMscModel();
142  msc3->SetHighEnergyLimit(highEnergyLimit);
143  msc3->SetRangeFactor(fRangeFactor);
144  msc3->SetGeomFactor(fGeomFactor);
145  msc3->SetSafetyFactor(fSafetyFactor);
146  msc3->SetLambdaLimit(fLambdaLimit);
147  msc3->SetStepLimitType(fStepLimitType);
148  msc3->SetLocked(true);
149  }
150 
151 #if G4VERSION_NUMBER >= 1110
152  G4TransportationWithMscType transportationWithMsc = G4EmParameters::Instance()->TransportationWithMsc();
153  if (transportationWithMsc != G4TransportationWithMscType::fDisabled) {
154  G4ProcessManager* procManager = particle->GetProcessManager();
155  // Remove default G4Transportation and replace with G4TransportationWithMsc.
156  G4VProcess* removed = procManager->RemoveProcess(0);
157  if (removed->GetProcessName() != "Transportation") {
158  G4Exception("CMSEmStandardPhysics::ConstructProcess",
159  "em0050",
160  FatalException,
161  "replaced process is not G4Transportation!");
162  }
163  G4TransportationWithMsc* transportWithMsc =
165  if (transportationWithMsc == G4TransportationWithMscType::fMultipleSteps) {
166  transportWithMsc->SetMultipleSteps(true);
167  }
168  transportWithMsc->AddMscModel(msc1);
169  transportWithMsc->AddMscModel(msc2);
170  if (nullptr != aRegion) {
171  transportWithMsc->AddMscModel(msc3, -1, aRegion);
172  }
173  if (nullptr != bRegion) {
174  transportWithMsc->AddMscModel(msc3, -1, bRegion);
175  }
176  procManager->AddProcess(transportWithMsc, -1, 0, 0);
177  } else
178 #endif
179  {
180  // Register as a separate process.
181  G4eMultipleScattering* msc = new G4eMultipleScattering;
182  msc->SetEmModel(msc1);
183  msc->SetEmModel(msc2);
184  if (nullptr != aRegion) {
185  msc->AddEmModel(-1, msc3, aRegion);
186  }
187  if (nullptr != bRegion) {
188  msc->AddEmModel(-1, msc3, bRegion);
189  }
190  ph->RegisterProcess(msc, particle);
191  }
192 
193  // single scattering
194  G4eCoulombScatteringModel* ssm = new G4eCoulombScatteringModel();
195  G4CoulombScattering* ss = new G4CoulombScattering();
196  ss->SetEmModel(ssm);
197  ss->SetMinKinEnergy(highEnergyLimit);
198  ssm->SetLowEnergyLimit(highEnergyLimit);
199  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
200 
201  ph->RegisterProcess(eioni, particle);
202  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
203  ph->RegisterProcess(ss, particle);
204 
205  // e+
206  particle = G4Positron::Positron();
207  eioni = new G4eIonisation();
208 
209  msc1 = new G4UrbanMscModel();
210  msc2 = new G4WentzelVIModel();
211  msc1->SetHighEnergyLimit(highEnergyLimit);
212  msc2->SetLowEnergyLimit(highEnergyLimit);
213 
214  // e-/e+ msc for HCAL and HGCAL using the Urban model
215  if (nullptr != aRegion || nullptr != bRegion) {
216  msc3 = new G4UrbanMscModel();
217  msc3->SetHighEnergyLimit(highEnergyLimit);
218  msc3->SetRangeFactor(fRangeFactor);
219  msc3->SetGeomFactor(fGeomFactor);
220  msc3->SetSafetyFactor(fSafetyFactor);
221  msc3->SetLambdaLimit(fLambdaLimit);
222  msc3->SetStepLimitType(fStepLimitType);
223  msc3->SetLocked(true);
224  }
225 
226 #if G4VERSION_NUMBER >= 1110
227  if (transportationWithMsc != G4TransportationWithMscType::fDisabled) {
228  G4ProcessManager* procManager = particle->GetProcessManager();
229  // Remove default G4Transportation and replace with G4TransportationWithMsc.
230  G4VProcess* removed = procManager->RemoveProcess(0);
231  if (removed->GetProcessName() != "Transportation") {
232  G4Exception("CMSEmStandardPhysics::ConstructProcess",
233  "em0050",
234  FatalException,
235  "replaced process is not G4Transportation!");
236  }
237  G4TransportationWithMsc* transportWithMsc =
239  if (transportationWithMsc == G4TransportationWithMscType::fMultipleSteps) {
240  transportWithMsc->SetMultipleSteps(true);
241  }
242  transportWithMsc->AddMscModel(msc1);
243  transportWithMsc->AddMscModel(msc2);
244  if (nullptr != aRegion) {
245  transportWithMsc->AddMscModel(msc3, -1, aRegion);
246  }
247  if (nullptr != bRegion) {
248  transportWithMsc->AddMscModel(msc3, -1, bRegion);
249  }
250  procManager->AddProcess(transportWithMsc, -1, 0, 0);
251  } else
252 #endif
253  {
254  // Register as a separate process.
255  G4eMultipleScattering* msc = new G4eMultipleScattering;
256  msc->SetEmModel(msc1);
257  msc->SetEmModel(msc2);
258  if (nullptr != aRegion) {
259  msc->AddEmModel(-1, msc3, aRegion);
260  }
261  if (nullptr != bRegion) {
262  msc->AddEmModel(-1, msc3, bRegion);
263  }
264  ph->RegisterProcess(msc, particle);
265  }
266 
267  // single scattering
268  ssm = new G4eCoulombScatteringModel();
269  ss = new G4CoulombScattering();
270  ss->SetEmModel(ssm);
271  ss->SetMinKinEnergy(highEnergyLimit);
272  ssm->SetLowEnergyLimit(highEnergyLimit);
273  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
274 
275  ph->RegisterProcess(eioni, particle);
276  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
277  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
278  ph->RegisterProcess(ss, particle);
279 
280  // generic ion
281  particle = G4GenericIon::GenericIon();
282  G4ionIonisation* ionIoni = new G4ionIonisation();
283  ph->RegisterProcess(hmsc, particle);
284  ph->RegisterProcess(ionIoni, particle);
285 
286  // muons, hadrons ions
287  G4EmBuilder::ConstructCharged(hmsc, pnuc);
288 }
Log< level::Info, true > LogVerbatim
void ConstructParticle() override
void ConstructProcess() override
G4MscStepLimitType fStepLimitType
CMSEmStandardPhysics(G4int ver, const edm::ParameterSet &p)