CMS 3D CMS Logo

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