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 47 of file CMSEmStandardPhysics.cc.

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

48  : G4VPhysicsConstructor("CMSEmStandard_emm") {
49  SetVerboseLevel(ver);
50  // EM parameters specific for this EM physics configuration
51  G4EmParameters* param = G4EmParameters::Instance();
52  param->SetDefaults();
53  param->SetVerbose(ver);
54  param->SetApplyCuts(true);
55  param->SetStepFunction(0.8, 1 * CLHEP::mm);
56  param->SetMscRangeFactor(0.2);
57  param->SetMscStepLimitType(fMinimal);
58  param->SetFluo(false);
59  SetPhysicsType(bElectromagnetic);
60  fRangeFactor = p.getParameter<double>("G4MscRangeFactor");
61  fGeomFactor = p.getParameter<double>("G4MscGeomFactor");
62  fSafetyFactor = p.getParameter<double>("G4MscSafetyFactor");
63  fLambdaLimit = p.getParameter<double>("G4MscLambdaLimit") * CLHEP::mm;
64  std::string msc = p.getParameter<std::string>("G4MscStepLimit");
65  fStepLimitType = fUseSafety;
66  if (msc == "UseSafetyPlus") {
67  fStepLimitType = fUseSafetyPlus;
68  }
69  if (msc == "Minimal") {
70  fStepLimitType = fMinimal;
71  }
72  double tcut = p.getParameter<double>("G4TrackingCut") * CLHEP::MeV;
73  param->SetLowestElectronEnergy(tcut);
74  param->SetLowestMuHadEnergy(tcut);
75  fG4HepEmActive = p.getParameter<bool>("G4HepEmActive");
76  if (fG4HepEmActive) {
77  // At the moment, G4HepEm supports only one configuration of MSC, so use
78  // the most generic parameters everywhere.
79  param->SetMscRangeFactor(fRangeFactor);
80  param->SetMscGeomFactor(fGeomFactor);
81  param->SetMscSafetyFactor(fSafetyFactor);
82  param->SetMscLambdaLimit(fLambdaLimit);
83  param->SetMscStepLimitType(fStepLimitType);
84  }
85 }
G4MscStepLimitType fStepLimitType

◆ ~CMSEmStandardPhysics()

CMSEmStandardPhysics::~CMSEmStandardPhysics ( )
overridedefault

Member Function Documentation

◆ ConstructParticle()

void CMSEmStandardPhysics::ConstructParticle ( )
override

Definition at line 87 of file CMSEmStandardPhysics.cc.

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

◆ ConstructProcess()

void CMSEmStandardPhysics::ConstructProcess ( )
override

Definition at line 92 of file CMSEmStandardPhysics.cc.

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

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