CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
CMSEmStandardPhysicsXS Class Reference

#include <CMSEmStandardPhysicsXS.h>

Inheritance diagram for CMSEmStandardPhysicsXS:

Public Member Functions

 CMSEmStandardPhysicsXS (G4int ver, const edm::ParameterSet &p)
 
void ConstructParticle () override
 
void ConstructProcess () override
 
 ~CMSEmStandardPhysicsXS () override
 

Private Attributes

G4double fGeomFactor
 
G4double fLambdaLimit
 
G4double fRangeFactor
 
G4double fSafetyFactor
 
G4MscStepLimitType fStepLimitType
 

Detailed Description

Definition at line 19 of file CMSEmStandardPhysicsXS.h.

Constructor & Destructor Documentation

◆ CMSEmStandardPhysicsXS()

CMSEmStandardPhysicsXS::CMSEmStandardPhysicsXS ( G4int  ver,
const edm::ParameterSet p 
)

Definition at line 63 of file CMSEmStandardPhysicsXS.cc.

References fGeomFactor, fLambdaLimit, fRangeFactor, fSafetyFactor, fStepLimitType, AlCaHLTBitMon_ParallelJobs::p, and AlCaHLTBitMon_QueryRunRegistry::string.

64  : G4VPhysicsConstructor("CMSEmStandard_emn") {
65  SetVerboseLevel(ver);
66  G4EmParameters* param = G4EmParameters::Instance();
67  param->SetDefaults();
68  param->SetVerbose(ver);
69  param->SetApplyCuts(true);
70  param->SetMinEnergy(100 * CLHEP::eV);
71  param->SetNumberOfBinsPerDecade(20);
72  param->SetStepFunction(0.8, 1 * CLHEP::mm);
73  param->SetMscRangeFactor(0.2);
74  param->SetMscStepLimitType(fMinimal);
75  param->SetFluo(true);
76  param->SetUseMottCorrection(true); // use Mott-correction for e-/e+ msc gs
77  SetPhysicsType(bElectromagnetic);
78  fRangeFactor = p.getParameter<double>("G4MscRangeFactor");
79  fGeomFactor = p.getParameter<double>("G4MscGeomFactor");
80  fSafetyFactor = p.getParameter<double>("G4MscSafetyFactor");
81  fLambdaLimit = p.getParameter<double>("G4MscLambdaLimit") * CLHEP::mm;
82  std::string msc = p.getParameter<std::string>("G4MscStepLimit");
83  fStepLimitType = fUseSafety;
84  if (msc == "UseSafetyPlus") {
85  fStepLimitType = fUseSafetyPlus;
86  }
87  if (msc == "Minimal") {
88  fStepLimitType = fMinimal;
89  }
90  double tcut = p.getParameter<double>("G4TrackingCut") * CLHEP::MeV;
91  param->SetLowestElectronEnergy(tcut);
92  param->SetLowestMuHadEnergy(tcut);
93 }
G4MscStepLimitType fStepLimitType

◆ ~CMSEmStandardPhysicsXS()

CMSEmStandardPhysicsXS::~CMSEmStandardPhysicsXS ( )
override

Definition at line 95 of file CMSEmStandardPhysicsXS.cc.

95 {}

Member Function Documentation

◆ ConstructParticle()

void CMSEmStandardPhysicsXS::ConstructParticle ( )
override

Definition at line 97 of file CMSEmStandardPhysicsXS.cc.

97  {
98  // minimal set of particles for EM physics
99  G4EmBuilder::ConstructMinimalEmSet();
100 }

◆ ConstructProcess()

void CMSEmStandardPhysicsXS::ConstructProcess ( )
override

Definition at line 102 of file CMSEmStandardPhysicsXS.cc.

References conv, callgraph::cs, nanoDQM_cfi::Electron, fGeomFactor, fLambdaLimit, fRangeFactor, fSafetyFactor, fStepLimitType, g4SimHits_cfi::G4GammaGeneralProcess, MaterialEffects_cfi::MultipleScattering, contentValuesCheck::ss, and HLT_2022v15_cff::verboseLevel.

102  {
103  if (verboseLevel > 0) {
104  edm::LogVerbatim("PhysicsList") << "### " << GetPhysicsName() << " Construct Processes";
105  }
106 
107  // This EM builder takes default models of Geant4 10 EMV.
108  // Multiple scattering by Urban for all particles
109  // except e+e- below 100 MeV for which the Urban model is used
110  G4EmBuilder::PrepareEMPhysics();
111 
112  G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
113 
114  // processes used by several particles
115  G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
116  G4NuclearStopping* pnuc(nullptr);
117 
118  // high energy limit for e+- scattering models
119  G4double highEnergyLimit = G4EmParameters::Instance()->MscEnergyLimit();
120  const G4Region* aRegion = G4RegionStore::GetInstance()->GetRegion("HcalRegion", false);
121  const G4Region* bRegion = G4RegionStore::GetInstance()->GetRegion("HGCalRegion", false);
122 
123  // Add gamma EM Processes
124  G4ParticleDefinition* particle = G4Gamma::Gamma();
125 
126  // Photoelectric
127  G4PhotoElectricEffect* pe = new G4PhotoElectricEffect();
128  G4VEmModel* theLivermorePEModel = new G4LivermorePhotoElectricModel();
129  pe->SetEmModel(theLivermorePEModel);
130 
131  // Compton scattering
132  G4ComptonScattering* cs = new G4ComptonScattering;
133  cs->SetEmModel(new G4KleinNishinaModel());
134 
135  // Gamma conversion
136  G4GammaConversion* gc = new G4GammaConversion();
137  G4VEmModel* conv = new G4BetheHeitler5DModel();
138  gc->SetEmModel(conv);
139 
140  if (G4EmParameters::Instance()->GeneralProcessActive()) {
142  sp->AddEmProcess(pe);
143  sp->AddEmProcess(cs);
144  sp->AddEmProcess(gc);
145  sp->AddEmProcess(new G4RayleighScattering());
146  G4LossTableManager::Instance()->SetGammaGeneralProcess(sp);
147  ph->RegisterProcess(sp, particle);
148  } else {
149  ph->RegisterProcess(pe, particle);
150  ph->RegisterProcess(cs, particle);
151  ph->RegisterProcess(gc, particle);
152  ph->RegisterProcess(new G4RayleighScattering(), particle);
153  }
154 
155  // e-
156  particle = G4Electron::Electron();
157 
158  // multiple scattering
159  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
160  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
161  msc1->SetHighEnergyLimit(highEnergyLimit);
162  msc2->SetLowEnergyLimit(highEnergyLimit);
163 
164  // msc for HCAL using the Urban model
165  G4UrbanMscModel* msc4 = nullptr;
166  if (nullptr != aRegion) {
167  msc4 = new G4UrbanMscModel();
168  msc4->SetHighEnergyLimit(highEnergyLimit);
169  msc4->SetRangeFactor(fRangeFactor);
170  msc4->SetGeomFactor(fGeomFactor);
171  msc4->SetSafetyFactor(fSafetyFactor);
172  msc4->SetLambdaLimit(fLambdaLimit);
173  msc4->SetStepLimitType(fStepLimitType);
174  msc4->SetLocked(true);
175  }
176 
177  // msc GS with Mott-correction
178  G4GoudsmitSaundersonMscModel* msc3 = nullptr;
179  if (nullptr != bRegion) {
180  msc3 = new G4GoudsmitSaundersonMscModel();
181  msc3->SetHighEnergyLimit(highEnergyLimit);
182  msc3->SetRangeFactor(0.08);
183  msc3->SetSkin(3.);
184  msc3->SetStepLimitType(fUseSafetyPlus);
185  msc3->SetLocked(true);
186  }
187 
188 #if G4VERSION_NUMBER >= 1110
189  G4TransportationWithMscType transportationWithMsc = G4EmParameters::Instance()->TransportationWithMsc();
190  if (transportationWithMsc != G4TransportationWithMscType::fDisabled) {
191  G4ProcessManager* procManager = particle->GetProcessManager();
192  // Remove default G4Transportation and replace with G4TransportationWithMsc.
193  G4VProcess* removed = procManager->RemoveProcess(0);
194  if (removed->GetProcessName() != "Transportation") {
195  G4Exception("CMSEmStandardPhysics::ConstructProcess",
196  "em0050",
197  FatalException,
198  "replaced process is not G4Transportation!");
199  }
200  G4TransportationWithMsc* transportWithMsc =
202  if (transportationWithMsc == G4TransportationWithMscType::fMultipleSteps) {
203  transportWithMsc->SetMultipleSteps(true);
204  }
205  transportWithMsc->AddMscModel(msc1);
206  transportWithMsc->AddMscModel(msc2);
207  if (msc4 != nullptr) {
208  transportWithMsc->AddMscModel(msc4, -1, aRegion);
209  }
210  if (msc3 != nullptr) {
211  transportWithMsc->AddMscModel(msc3, -1, bRegion);
212  }
213  procManager->AddProcess(transportWithMsc, -1, 0, 0);
214  } else
215 #endif
216  {
217  // Register as a separate process.
218  G4eMultipleScattering* msc = new G4eMultipleScattering;
219  msc->SetEmModel(msc1);
220  msc->SetEmModel(msc2);
221  if (msc4 != nullptr) {
222  msc->AddEmModel(-1, msc4, aRegion);
223  }
224  if (msc3 != nullptr) {
225  msc->AddEmModel(-1, msc3, bRegion);
226  }
227  ph->RegisterProcess(msc, particle);
228  }
229 
230  // single scattering
231  G4eCoulombScatteringModel* ssm = new G4eCoulombScatteringModel();
232  G4CoulombScattering* ss = new G4CoulombScattering();
233  ss->SetEmModel(ssm);
234  ss->SetMinKinEnergy(highEnergyLimit);
235  ssm->SetLowEnergyLimit(highEnergyLimit);
236  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
237 
238  // ionisation
239  G4eIonisation* eioni = new G4eIonisation();
240 
241  // bremsstrahlung
242  G4eBremsstrahlung* brem = new G4eBremsstrahlung();
243  G4SeltzerBergerModel* br1 = new G4SeltzerBergerModel();
244  G4eBremsstrahlungRelModel* br2 = new G4eBremsstrahlungRelModel();
245  br1->SetAngularDistribution(new G4Generator2BS());
246  br2->SetAngularDistribution(new G4Generator2BS());
247  brem->SetEmModel(br1);
248  brem->SetEmModel(br2);
249  br1->SetHighEnergyLimit(CLHEP::GeV);
250 
251  G4ePairProduction* ee = new G4ePairProduction();
252 
253  // register processes
254  ph->RegisterProcess(eioni, particle);
255  ph->RegisterProcess(brem, particle);
256  ph->RegisterProcess(ee, particle);
257  ph->RegisterProcess(ss, particle);
258 
259  // e+
260  particle = G4Positron::Positron();
261 
262  // multiple scattering
263  msc1 = new G4UrbanMscModel();
264  msc2 = new G4WentzelVIModel();
265  msc1->SetHighEnergyLimit(highEnergyLimit);
266  msc2->SetLowEnergyLimit(highEnergyLimit);
267 
268  // msc for HCAL using the Urban model
269  if (nullptr != aRegion) {
270  msc4 = new G4UrbanMscModel();
271  msc4->SetHighEnergyLimit(highEnergyLimit);
272  msc4->SetRangeFactor(fRangeFactor);
273  msc4->SetGeomFactor(fGeomFactor);
274  msc4->SetSafetyFactor(fSafetyFactor);
275  msc4->SetLambdaLimit(fLambdaLimit);
276  msc4->SetStepLimitType(fStepLimitType);
277  msc4->SetLocked(true);
278  }
279 
280  // msc GS with Mott-correction
281  if (nullptr != bRegion) {
282  msc3 = new G4GoudsmitSaundersonMscModel();
283  msc3->SetHighEnergyLimit(highEnergyLimit);
284  msc3->SetRangeFactor(0.08);
285  msc3->SetSkin(3.);
286  msc3->SetStepLimitType(fUseSafetyPlus);
287  msc3->SetLocked(true);
288  }
289 
290 #if G4VERSION_NUMBER >= 1110
291  if (transportationWithMsc != G4TransportationWithMscType::fDisabled) {
292  G4ProcessManager* procManager = particle->GetProcessManager();
293  // Remove default G4Transportation and replace with G4TransportationWithMsc.
294  G4VProcess* removed = procManager->RemoveProcess(0);
295  if (removed->GetProcessName() != "Transportation") {
296  G4Exception("CMSEmStandardPhysics::ConstructProcess",
297  "em0050",
298  FatalException,
299  "replaced process is not G4Transportation!");
300  }
301  G4TransportationWithMsc* transportWithMsc =
303  if (transportationWithMsc == G4TransportationWithMscType::fMultipleSteps) {
304  transportWithMsc->SetMultipleSteps(true);
305  }
306  transportWithMsc->AddMscModel(msc1);
307  transportWithMsc->AddMscModel(msc2);
308  if (msc4 != nullptr) {
309  transportWithMsc->AddMscModel(msc4, -1, aRegion);
310  }
311  if (msc3 != nullptr) {
312  transportWithMsc->AddMscModel(msc3, -1, bRegion);
313  }
314  procManager->AddProcess(transportWithMsc, -1, 0, 0);
315  } else
316 #endif
317  {
318  // Register as a separate process.
319  G4eMultipleScattering* msc = new G4eMultipleScattering;
320  msc->SetEmModel(msc1);
321  msc->SetEmModel(msc2);
322  if (msc4 != nullptr) {
323  msc->AddEmModel(-1, msc4, aRegion);
324  }
325  if (msc3 != nullptr) {
326  msc->AddEmModel(-1, msc3, bRegion);
327  }
328  ph->RegisterProcess(msc, particle);
329  }
330 
331  // single scattering
332  ssm = new G4eCoulombScatteringModel();
333  ss = new G4CoulombScattering();
334  ss->SetEmModel(ssm);
335  ss->SetMinKinEnergy(highEnergyLimit);
336  ssm->SetLowEnergyLimit(highEnergyLimit);
337  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
338 
339  // ionisation
340  eioni = new G4eIonisation();
341 
342  // bremsstrahlung
343  brem = new G4eBremsstrahlung();
344  br1 = new G4SeltzerBergerModel();
345  br2 = new G4eBremsstrahlungRelModel();
346  br1->SetAngularDistribution(new G4Generator2BS());
347  br2->SetAngularDistribution(new G4Generator2BS());
348  brem->SetEmModel(br1);
349  brem->SetEmModel(br2);
350  br1->SetHighEnergyLimit(CLHEP::GeV);
351 
352  // register processes
353  ph->RegisterProcess(eioni, particle);
354  ph->RegisterProcess(brem, particle);
355  ph->RegisterProcess(ee, particle);
356  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
357  ph->RegisterProcess(ss, particle);
358 
359  // generic ion
360  particle = G4GenericIon::GenericIon();
361  G4ionIonisation* ionIoni = new G4ionIonisation();
362  ph->RegisterProcess(hmsc, particle);
363  ph->RegisterProcess(ionIoni, particle);
364 
365  // muons, hadrons, ions
366  G4EmBuilder::ConstructCharged(hmsc, pnuc);
367 }
Log< level::Info, true > LogVerbatim
EPOS::IO_EPOS conv
G4MscStepLimitType fStepLimitType

Member Data Documentation

◆ fGeomFactor

G4double CMSEmStandardPhysicsXS::fGeomFactor
private

Definition at line 29 of file CMSEmStandardPhysicsXS.h.

Referenced by CMSEmStandardPhysicsXS(), and ConstructProcess().

◆ fLambdaLimit

G4double CMSEmStandardPhysicsXS::fLambdaLimit
private

Definition at line 31 of file CMSEmStandardPhysicsXS.h.

Referenced by CMSEmStandardPhysicsXS(), and ConstructProcess().

◆ fRangeFactor

G4double CMSEmStandardPhysicsXS::fRangeFactor
private

Definition at line 28 of file CMSEmStandardPhysicsXS.h.

Referenced by CMSEmStandardPhysicsXS(), and ConstructProcess().

◆ fSafetyFactor

G4double CMSEmStandardPhysicsXS::fSafetyFactor
private

Definition at line 30 of file CMSEmStandardPhysicsXS.h.

Referenced by CMSEmStandardPhysicsXS(), and ConstructProcess().

◆ fStepLimitType

G4MscStepLimitType CMSEmStandardPhysicsXS::fStepLimitType
private

Definition at line 32 of file CMSEmStandardPhysicsXS.h.

Referenced by CMSEmStandardPhysicsXS(), and ConstructProcess().