CMS 3D CMS Logo

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

#include <CMSEmStandardPhysics.h>

Inheritance diagram for CMSEmStandardPhysics:

Public Member Functions

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

Private Attributes

bool fG4HepEmActive
 
G4double fGeomFactor
 
G4double fLambdaLimit
 
G4double fRangeFactor
 
G4double fSafetyFactor
 
G4MscStepLimitType fStepLimitType
 

Detailed Description

Definition at line 17 of file CMSEmStandardPhysics.h.

Constructor & Destructor Documentation

◆ CMSEmStandardPhysics()

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

Definition at line 49 of file CMSEmStandardPhysics.cc.

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

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

◆ ~CMSEmStandardPhysics()

CMSEmStandardPhysics::~CMSEmStandardPhysics ( )
overridedefault

Member Function Documentation

◆ ConstructParticle()

void CMSEmStandardPhysics::ConstructParticle ( )
override

Definition at line 89 of file CMSEmStandardPhysics.cc.

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

◆ ConstructProcess()

void CMSEmStandardPhysics::ConstructProcess ( )
override

Definition at line 94 of file CMSEmStandardPhysics.cc.

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

94  {
95  if (verboseLevel > 0) {
96  edm::LogVerbatim("PhysicsList") << "### " << GetPhysicsName() << " Construct EM Processes";
97  }
98 
99  // This EM builder takes default models of Geant4 10 EMV.
100  // Multiple scattering by WentzelVI for all particles except:
101  // a) e+e- below 100 MeV for which the Urban model is used
102  // b) ions for which Urban model is used
103  G4EmBuilder::PrepareEMPhysics();
104 
105  G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
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  G4eIonisation* eioni = new G4eIonisation();
140 
141  G4UrbanMscModel* msc1 = new G4UrbanMscModel();
142  G4WentzelVIModel* msc2 = new G4WentzelVIModel();
143  msc1->SetHighEnergyLimit(highEnergyLimit);
144  msc2->SetLowEnergyLimit(highEnergyLimit);
145 
146  // e-/e+ msc for HCAL and HGCAL using the Urban model
147  G4UrbanMscModel* msc3 = nullptr;
148  if (nullptr != aRegion || nullptr != bRegion) {
149  msc3 = new G4UrbanMscModel();
150  msc3->SetHighEnergyLimit(highEnergyLimit);
151  msc3->SetRangeFactor(fRangeFactor);
152  msc3->SetGeomFactor(fGeomFactor);
153  msc3->SetSafetyFactor(fSafetyFactor);
154  msc3->SetLambdaLimit(fLambdaLimit);
155  msc3->SetStepLimitType(fStepLimitType);
156  msc3->SetLocked(true);
157  }
158 
159  G4TransportationWithMscType transportationWithMsc = param->TransportationWithMsc();
160  if (transportationWithMsc != G4TransportationWithMscType::fDisabled) {
161  // Remove default G4Transportation and replace with G4TransportationWithMsc.
162  G4ProcessManager* procManager = particle->GetProcessManager();
163  G4VProcess* removed = procManager->RemoveProcess(0);
164  if (removed->GetProcessName() != "Transportation") {
165  G4Exception("CMSEmStandardPhysics::ConstructProcess",
166  "em0050",
167  FatalException,
168  "replaced process is not G4Transportation!");
169  }
170  G4TransportationWithMsc* transportWithMsc =
172  if (transportationWithMsc == G4TransportationWithMscType::fMultipleSteps) {
173  transportWithMsc->SetMultipleSteps(true);
174  }
175  transportWithMsc->AddMscModel(msc1);
176  transportWithMsc->AddMscModel(msc2);
177  if (nullptr != aRegion) {
178  transportWithMsc->AddMscModel(msc3, -1, aRegion);
179  }
180  if (nullptr != bRegion) {
181  transportWithMsc->AddMscModel(msc3, -1, bRegion);
182  }
183  procManager->AddProcess(transportWithMsc, -1, 0, 0);
184  } else {
185  // Multiple scattering is registered as a separate process
186  G4eMultipleScattering* msc = new G4eMultipleScattering;
187  msc->SetEmModel(msc1);
188  msc->SetEmModel(msc2);
189  if (nullptr != aRegion) {
190  msc->AddEmModel(-1, msc3, aRegion);
191  }
192  if (nullptr != bRegion) {
193  msc->AddEmModel(-1, msc3, bRegion);
194  }
195  ph->RegisterProcess(msc, particle);
196  }
197 
198  // single scattering
199  G4eCoulombScatteringModel* ssm = new G4eCoulombScatteringModel();
200  G4CoulombScattering* ss = new G4CoulombScattering();
201  ss->SetEmModel(ssm);
202  ss->SetMinKinEnergy(highEnergyLimit);
203  ssm->SetLowEnergyLimit(highEnergyLimit);
204  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
205 
206  ph->RegisterProcess(eioni, particle);
207  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
208  ph->RegisterProcess(ss, particle);
209 
210  // e+
211  particle = G4Positron::Positron();
212  eioni = new G4eIonisation();
213 
214  msc1 = new G4UrbanMscModel();
215  msc2 = new G4WentzelVIModel();
216  msc1->SetHighEnergyLimit(highEnergyLimit);
217  msc2->SetLowEnergyLimit(highEnergyLimit);
218 
219  // e-/e+ msc for HCAL and HGCAL using the Urban model
220  if (nullptr != aRegion || nullptr != bRegion) {
221  msc3 = new G4UrbanMscModel();
222  msc3->SetHighEnergyLimit(highEnergyLimit);
223  msc3->SetRangeFactor(fRangeFactor);
224  msc3->SetGeomFactor(fGeomFactor);
225  msc3->SetSafetyFactor(fSafetyFactor);
226  msc3->SetLambdaLimit(fLambdaLimit);
227  msc3->SetStepLimitType(fStepLimitType);
228  msc3->SetLocked(true);
229  }
230 
231  if (transportationWithMsc != G4TransportationWithMscType::fDisabled) {
232  G4ProcessManager* procManager = particle->GetProcessManager();
233  // Remove default G4Transportation and replace with G4TransportationWithMsc.
234  G4VProcess* removed = procManager->RemoveProcess(0);
235  if (removed->GetProcessName() != "Transportation") {
236  G4Exception("CMSEmStandardPhysics::ConstructProcess",
237  "em0050",
238  FatalException,
239  "replaced process is not G4Transportation!");
240  }
241  G4TransportationWithMsc* transportWithMsc =
243  if (transportationWithMsc == G4TransportationWithMscType::fMultipleSteps) {
244  transportWithMsc->SetMultipleSteps(true);
245  }
246  transportWithMsc->AddMscModel(msc1);
247  transportWithMsc->AddMscModel(msc2);
248  if (nullptr != aRegion) {
249  transportWithMsc->AddMscModel(msc3, -1, aRegion);
250  }
251  if (nullptr != bRegion) {
252  transportWithMsc->AddMscModel(msc3, -1, bRegion);
253  }
254  procManager->AddProcess(transportWithMsc, -1, 0, 0);
255  } else {
256  // Register as a separate process.
257  G4eMultipleScattering* msc = new G4eMultipleScattering;
258  msc->SetEmModel(msc1);
259  msc->SetEmModel(msc2);
260  if (nullptr != aRegion) {
261  msc->AddEmModel(-1, msc3, aRegion);
262  }
263  if (nullptr != bRegion) {
264  msc->AddEmModel(-1, msc3, bRegion);
265  }
266  ph->RegisterProcess(msc, particle);
267  }
268 
269  // single scattering
270  ssm = new G4eCoulombScatteringModel();
271  ss = new G4CoulombScattering();
272  ss->SetEmModel(ssm);
273  ss->SetMinKinEnergy(highEnergyLimit);
274  ssm->SetLowEnergyLimit(highEnergyLimit);
275  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
276 
277  ph->RegisterProcess(eioni, particle);
278  ph->RegisterProcess(new G4eBremsstrahlung(), particle);
279  ph->RegisterProcess(new G4eplusAnnihilation(), particle);
280  ph->RegisterProcess(ss, particle);
281 
282  if (fG4HepEmActive) {
283  auto* hepEmTM = new CMSHepEmTrackingManager(highEnergyLimit);
284  G4Electron::Electron()->SetTrackingManager(hepEmTM);
285  G4Positron::Positron()->SetTrackingManager(hepEmTM);
286  }
287 
288  // generic ion
289  particle = G4GenericIon::GenericIon();
290  G4ionIonisation* ionIoni = new G4ionIonisation();
291  ph->RegisterProcess(hmsc, particle);
292  ph->RegisterProcess(ionIoni, particle);
293 
294  // muons, hadrons ions
295  G4EmBuilder::ConstructCharged(hmsc, pnuc);
296 }
Log< level::Info, true > LogVerbatim
G4MscStepLimitType fStepLimitType

Member Data Documentation

◆ fG4HepEmActive

bool CMSEmStandardPhysics::fG4HepEmActive
private

Definition at line 31 of file CMSEmStandardPhysics.h.

Referenced by CMSEmStandardPhysics(), and ConstructProcess().

◆ fGeomFactor

G4double CMSEmStandardPhysics::fGeomFactor
private

Definition at line 27 of file CMSEmStandardPhysics.h.

Referenced by CMSEmStandardPhysics(), and ConstructProcess().

◆ fLambdaLimit

G4double CMSEmStandardPhysics::fLambdaLimit
private

Definition at line 29 of file CMSEmStandardPhysics.h.

Referenced by CMSEmStandardPhysics(), and ConstructProcess().

◆ fRangeFactor

G4double CMSEmStandardPhysics::fRangeFactor
private

Definition at line 26 of file CMSEmStandardPhysics.h.

Referenced by CMSEmStandardPhysics(), and ConstructProcess().

◆ fSafetyFactor

G4double CMSEmStandardPhysics::fSafetyFactor
private

Definition at line 28 of file CMSEmStandardPhysics.h.

Referenced by CMSEmStandardPhysics(), and ConstructProcess().

◆ fStepLimitType

G4MscStepLimitType CMSEmStandardPhysics::fStepLimitType
private

Definition at line 30 of file CMSEmStandardPhysics.h.

Referenced by CMSEmStandardPhysics(), and ConstructProcess().