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=default
 

Private Attributes

bool fG4HepEmActive
 
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 48 of file CMSEmStandardPhysicsXS.cc.

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

49  : G4VPhysicsConstructor("CMSEmStandard_emn") {
50  SetVerboseLevel(ver);
51  // EM parameters specific for this EM physics configuration
52  G4EmParameters* param = G4EmParameters::Instance();
53  param->SetDefaults();
54  param->SetVerbose(ver);
55  param->SetApplyCuts(true);
56  param->SetStepFunction(0.8, 1 * CLHEP::mm);
57  param->SetMscRangeFactor(0.2);
58  param->SetMscStepLimitType(fMinimal);
59  param->SetFluo(false);
60  param->SetUseMottCorrection(true); // use Mott-correction for e-/e+ msc gs
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  if (fG4HepEmActive) {
79  // At the moment, G4HepEm supports only one configuration of MSC, so use
80  // the most generic parameters everywhere.
81  param->SetMscRangeFactor(fRangeFactor);
82  param->SetMscGeomFactor(fGeomFactor);
83  param->SetMscSafetyFactor(fSafetyFactor);
84  param->SetMscLambdaLimit(fLambdaLimit);
85  param->SetMscStepLimitType(fStepLimitType);
86  }
87 }
G4MscStepLimitType fStepLimitType

◆ ~CMSEmStandardPhysicsXS()

CMSEmStandardPhysicsXS::~CMSEmStandardPhysicsXS ( )
overridedefault

Member Function Documentation

◆ ConstructParticle()

void CMSEmStandardPhysicsXS::ConstructParticle ( )
override

Definition at line 89 of file CMSEmStandardPhysicsXS.cc.

89  {
90  // minimal set of particles for EM physics
91  G4EmBuilder::ConstructMinimalEmSet();
92 }

◆ ConstructProcess()

void CMSEmStandardPhysicsXS::ConstructProcess ( )
override

Definition at line 94 of file CMSEmStandardPhysicsXS.cc.

References nanoDQM_cfi::Electron, fG4HepEmActive, fGeomFactor, fLambdaLimit, fRangeFactor, fSafetyFactor, fStepLimitType, g4SimHits_cfi::G4GammaGeneralProcess, MaterialEffects_cfi::MultipleScattering, contentValuesCheck::ss, and HLT_2024v11_cff::verboseLevel.

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

Member Data Documentation

◆ fG4HepEmActive

bool CMSEmStandardPhysicsXS::fG4HepEmActive
private

Definition at line 33 of file CMSEmStandardPhysicsXS.h.

Referenced by CMSEmStandardPhysicsXS(), and ConstructProcess().

◆ 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().