4 #include "G4SystemOfUnits.hh" 5 #include "G4ParticleDefinition.hh" 6 #include "G4LossTableManager.hh" 7 #include "G4EmParameters.hh" 8 #include "G4EmBuilder.hh" 10 #include "G4ComptonScattering.hh" 11 #include "G4GammaConversion.hh" 12 #include "G4PhotoElectricEffect.hh" 13 #include "G4RayleighScattering.hh" 14 #include "G4PEEffectFluoModel.hh" 15 #include "G4KleinNishinaModel.hh" 16 #include "G4LowEPComptonModel.hh" 17 #include "G4BetheHeitler5DModel.hh" 18 #include "G4LivermorePhotoElectricModel.hh" 20 #include "G4eMultipleScattering.hh" 21 #include "G4hMultipleScattering.hh" 22 #include "G4MscStepLimitType.hh" 23 #include "G4UrbanMscModel.hh" 24 #include "G4GoudsmitSaundersonMscModel.hh" 25 #include "G4DummyModel.hh" 26 #include "G4WentzelVIModel.hh" 27 #include "G4CoulombScattering.hh" 28 #include "G4eCoulombScatteringModel.hh" 30 #include "G4eIonisation.hh" 31 #include "G4eBremsstrahlung.hh" 32 #include "G4Generator2BS.hh" 33 #include "G4SeltzerBergerModel.hh" 34 #include "G4ePairProduction.hh" 35 #include "G4UniversalFluctuation.hh" 37 #include "G4eplusAnnihilation.hh" 39 #include "G4hIonisation.hh" 40 #include "G4ionIonisation.hh" 43 #include "G4Electron.hh" 44 #include "G4Positron.hh" 45 #include "G4GenericIon.hh" 47 #include "G4PhysicsListHelper.hh" 48 #include "G4BuilderType.hh" 49 #include "G4GammaGeneralProcess.hh" 51 #include "G4RegionStore.hh" 52 #include "G4Region.hh" 53 #include "G4GammaGeneralProcess.hh" 55 #include "G4SystemOfUnits.hh" 58 : G4VPhysicsConstructor(
"CMSEmStandard_emn") {
60 G4EmParameters* param = G4EmParameters::Instance();
62 param->SetVerbose(ver);
63 param->SetApplyCuts(
true);
64 param->SetMinEnergy(100 * CLHEP::eV);
65 param->SetNumberOfBinsPerDecade(20);
66 param->SetStepFunction(0.8, 1 * CLHEP::mm);
67 param->SetMscRangeFactor(0.2);
68 param->SetMscStepLimitType(fMinimal);
70 param->SetUseMottCorrection(
true);
71 SetPhysicsType(bElectromagnetic);
75 fLambdaLimit =
p.getParameter<
double>(
"G4MscLambdaLimit") * CLHEP::mm;
78 if (msc ==
"UseSafetyPlus") {
81 if (msc ==
"Minimal") {
90 G4EmBuilder::ConstructMinimalEmSet();
95 edm::LogVerbatim(
"PhysicsList") <<
"### " << GetPhysicsName() <<
" Construct Processes";
101 G4EmBuilder::PrepareEMPhysics();
103 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
106 G4hMultipleScattering* hmsc =
new G4hMultipleScattering(
"ionmsc");
107 G4NuclearStopping* pnuc(
nullptr);
110 G4double highEnergyLimit = G4EmParameters::Instance()->MscEnergyLimit();
111 const G4Region* aRegion = G4RegionStore::GetInstance()->GetRegion(
"HcalRegion",
false);
112 const G4Region* bRegion = G4RegionStore::GetInstance()->GetRegion(
"HGCalRegion",
false);
115 G4ParticleDefinition* particle = G4Gamma::Gamma();
118 G4PhotoElectricEffect* pe =
new G4PhotoElectricEffect();
119 G4VEmModel* theLivermorePEModel =
new G4LivermorePhotoElectricModel();
120 pe->SetEmModel(theLivermorePEModel);
123 G4ComptonScattering*
cs =
new G4ComptonScattering;
124 cs->SetEmModel(
new G4KleinNishinaModel());
127 G4GammaConversion* gc =
new G4GammaConversion();
128 G4VEmModel*
conv =
new G4BetheHeitler5DModel();
129 gc->SetEmModel(
conv);
131 if (G4EmParameters::Instance()->GeneralProcessActive()) {
132 G4GammaGeneralProcess* sp =
new G4GammaGeneralProcess();
133 sp->AddEmProcess(pe);
134 sp->AddEmProcess(
cs);
135 sp->AddEmProcess(gc);
136 sp->AddEmProcess(
new G4RayleighScattering());
137 G4LossTableManager::Instance()->SetGammaGeneralProcess(sp);
138 ph->RegisterProcess(sp, particle);
140 ph->RegisterProcess(pe, particle);
141 ph->RegisterProcess(
cs, particle);
142 ph->RegisterProcess(gc, particle);
143 ph->RegisterProcess(
new G4RayleighScattering(), particle);
150 G4eMultipleScattering* msc =
new G4eMultipleScattering();
151 G4UrbanMscModel* msc1 =
new G4UrbanMscModel();
152 G4WentzelVIModel* msc2 =
new G4WentzelVIModel();
153 msc1->SetHighEnergyLimit(highEnergyLimit);
154 msc2->SetLowEnergyLimit(highEnergyLimit);
155 msc->SetEmModel(msc1);
156 msc->SetEmModel(msc2);
159 if (
nullptr != aRegion) {
160 G4UrbanMscModel* msc4 =
new G4UrbanMscModel();
161 msc4->SetHighEnergyLimit(highEnergyLimit);
167 msc4->SetLocked(
true);
168 msc->AddEmModel(-1, msc4, aRegion);
172 if (
nullptr != bRegion) {
173 G4GoudsmitSaundersonMscModel* msc3 =
new G4GoudsmitSaundersonMscModel();
174 msc3->SetHighEnergyLimit(highEnergyLimit);
175 msc3->SetRangeFactor(0.08);
177 msc3->SetStepLimitType(fUseSafetyPlus);
178 msc3->SetLocked(
true);
179 msc->AddEmModel(-1, msc3, bRegion);
183 G4eCoulombScatteringModel* ssm =
new G4eCoulombScatteringModel();
184 G4CoulombScattering*
ss =
new G4CoulombScattering();
186 ss->SetMinKinEnergy(highEnergyLimit);
187 ssm->SetLowEnergyLimit(highEnergyLimit);
188 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
191 G4eIonisation* eioni =
new G4eIonisation();
194 G4eBremsstrahlung* brem =
new G4eBremsstrahlung();
195 G4SeltzerBergerModel* br1 =
new G4SeltzerBergerModel();
196 G4eBremsstrahlungRelModel* br2 =
new G4eBremsstrahlungRelModel();
197 br1->SetAngularDistribution(
new G4Generator2BS());
198 br2->SetAngularDistribution(
new G4Generator2BS());
199 brem->SetEmModel(br1);
200 brem->SetEmModel(br2);
201 br1->SetHighEnergyLimit(CLHEP::GeV);
203 G4ePairProduction* ee =
new G4ePairProduction();
206 ph->RegisterProcess(msc, particle);
207 ph->RegisterProcess(eioni, particle);
208 ph->RegisterProcess(brem, particle);
209 ph->RegisterProcess(ee, particle);
210 ph->RegisterProcess(
ss, particle);
213 particle = G4Positron::Positron();
216 msc =
new G4eMultipleScattering();
217 msc1 =
new G4UrbanMscModel();
218 msc2 =
new G4WentzelVIModel();
219 msc1->SetHighEnergyLimit(highEnergyLimit);
220 msc2->SetLowEnergyLimit(highEnergyLimit);
221 msc->SetEmModel(msc1);
222 msc->SetEmModel(msc2);
225 if (
nullptr != aRegion) {
226 G4UrbanMscModel* msc4 =
new G4UrbanMscModel();
227 msc4->SetHighEnergyLimit(highEnergyLimit);
233 msc4->SetLocked(
true);
234 msc->AddEmModel(-1, msc4, aRegion);
238 if (
nullptr != bRegion) {
239 G4GoudsmitSaundersonMscModel* msc3 =
new G4GoudsmitSaundersonMscModel();
240 msc3->SetHighEnergyLimit(highEnergyLimit);
241 msc3->SetRangeFactor(0.08);
243 msc3->SetStepLimitType(fUseSafetyPlus);
244 msc3->SetLocked(
true);
245 msc->AddEmModel(-1, msc3, bRegion);
249 ssm =
new G4eCoulombScatteringModel();
250 ss =
new G4CoulombScattering();
252 ss->SetMinKinEnergy(highEnergyLimit);
253 ssm->SetLowEnergyLimit(highEnergyLimit);
254 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
257 eioni =
new G4eIonisation();
260 brem =
new G4eBremsstrahlung();
261 br1 =
new G4SeltzerBergerModel();
262 br2 =
new G4eBremsstrahlungRelModel();
263 br1->SetAngularDistribution(
new G4Generator2BS());
264 br2->SetAngularDistribution(
new G4Generator2BS());
265 brem->SetEmModel(br1);
266 brem->SetEmModel(br2);
267 br1->SetHighEnergyLimit(CLHEP::GeV);
270 ph->RegisterProcess(msc, particle);
271 ph->RegisterProcess(eioni, particle);
272 ph->RegisterProcess(brem, particle);
273 ph->RegisterProcess(ee, particle);
274 ph->RegisterProcess(
new G4eplusAnnihilation(), particle);
275 ph->RegisterProcess(
ss, particle);
278 particle = G4GenericIon::GenericIon();
279 G4ionIonisation* ionIoni =
new G4ionIonisation();
280 ph->RegisterProcess(hmsc, particle);
281 ph->RegisterProcess(ionIoni, particle);
284 G4EmBuilder::ConstructCharged(hmsc, pnuc);
Log< level::Info, true > LogVerbatim
void ConstructProcess() override
~CMSEmStandardPhysicsXS() override
CMSEmStandardPhysicsXS(G4int ver, const edm::ParameterSet &p)
void ConstructParticle() override
G4MscStepLimitType fStepLimitType