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 57 of file CMSEmStandardPhysicsXS.cc.

58  : G4VPhysicsConstructor("CMSEmStandard_emn") {
59  SetVerboseLevel(ver);
60  G4EmParameters* param = G4EmParameters::Instance();
61  param->SetDefaults();
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);
69  param->SetFluo(true);
70  param->SetUseMottCorrection(true); // use Mott-correction for e-/e+ msc gs
71  SetPhysicsType(bElectromagnetic);
72  fRangeFactor = p.getParameter<double>("G4MscRangeFactor");
73  fGeomFactor = p.getParameter<double>("G4MscGeomFactor");
74  fSafetyFactor = p.getParameter<double>("G4MscSafetyFactor");
75  fLambdaLimit = p.getParameter<double>("G4MscLambdaLimit") * CLHEP::mm;
76  std::string msc = p.getParameter<std::string>("G4MscStepLimit");
77  fStepLimitType = fUseSafety;
78  if (msc == "UseSafetyPlus") {
79  fStepLimitType = fUseSafetyPlus;
80  }
81  if (msc == "Minimal") {
82  fStepLimitType = fMinimal;
83  }
84 }

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

◆ ~CMSEmStandardPhysicsXS()

CMSEmStandardPhysicsXS::~CMSEmStandardPhysicsXS ( )
override

Definition at line 86 of file CMSEmStandardPhysicsXS.cc.

86 {}

Member Function Documentation

◆ ConstructParticle()

void CMSEmStandardPhysicsXS::ConstructParticle ( )
override

Definition at line 88 of file CMSEmStandardPhysicsXS.cc.

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

◆ ConstructProcess()

void CMSEmStandardPhysicsXS::ConstructProcess ( )
override

Definition at line 93 of file CMSEmStandardPhysicsXS.cc.

93  {
94  if (verboseLevel > 0) {
95  edm::LogVerbatim("PhysicsList") << "### " << GetPhysicsName() << " Construct Processes";
96  }
97 
98  // This EM builder takes default models of Geant4 10 EMV.
99  // Multiple scattering by Urban for all particles
100  // except e+e- below 100 MeV for which the Urban model is used
101  G4EmBuilder::PrepareEMPhysics();
102 
103  G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
104 
105  // processes used by several particles
106  G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
107  G4NuclearStopping* pnuc(nullptr);
108 
109  // high energy limit for e+- scattering models
110  G4double highEnergyLimit = G4EmParameters::Instance()->MscEnergyLimit();
111  const G4Region* aRegion = G4RegionStore::GetInstance()->GetRegion("HcalRegion", false);
112  const G4Region* bRegion = G4RegionStore::GetInstance()->GetRegion("HGCalRegion", false);
113 
114  // Add gamma EM Processes
115  G4ParticleDefinition* particle = G4Gamma::Gamma();
116 
117  // Photoelectric
118  G4PhotoElectricEffect* pe = new G4PhotoElectricEffect();
119  G4VEmModel* theLivermorePEModel = new G4LivermorePhotoElectricModel();
120  pe->SetEmModel(theLivermorePEModel);
121 
122  // Compton scattering
123  G4ComptonScattering* cs = new G4ComptonScattering;
124  cs->SetEmModel(new G4KleinNishinaModel());
125 
126  // Gamma conversion
127  G4GammaConversion* gc = new G4GammaConversion();
128  G4VEmModel* conv = new G4BetheHeitler5DModel();
129  gc->SetEmModel(conv);
130 
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);
139  } else {
140  ph->RegisterProcess(pe, particle);
141  ph->RegisterProcess(cs, particle);
142  ph->RegisterProcess(gc, particle);
143  ph->RegisterProcess(new G4RayleighScattering(), particle);
144  }
145 
146  // e-
147  particle = G4Electron::Electron();
148 
149  // multiple scattering
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);
157 
158  // msc for HCAL using the Urban model
159  if (nullptr != aRegion) {
160  G4UrbanMscModel* msc4 = new G4UrbanMscModel();
161  msc4->SetHighEnergyLimit(highEnergyLimit);
162  msc4->SetRangeFactor(fRangeFactor);
163  msc4->SetGeomFactor(fGeomFactor);
164  msc4->SetSafetyFactor(fSafetyFactor);
165  msc4->SetLambdaLimit(fLambdaLimit);
166  msc4->SetStepLimitType(fStepLimitType);
167  msc4->SetLocked(true);
168  msc->AddEmModel(-1, msc4, aRegion);
169  }
170 
171  // msc GS with Mott-correction
172  if (nullptr != bRegion) {
173  G4GoudsmitSaundersonMscModel* msc3 = new G4GoudsmitSaundersonMscModel();
174  msc3->SetHighEnergyLimit(highEnergyLimit);
175  msc3->SetRangeFactor(0.08);
176  msc3->SetSkin(3.);
177  msc3->SetStepLimitType(fUseSafetyPlus);
178  msc3->SetLocked(true);
179  msc->AddEmModel(-1, msc3, bRegion);
180  }
181 
182  // single scattering
183  G4eCoulombScatteringModel* ssm = new G4eCoulombScatteringModel();
184  G4CoulombScattering* ss = new G4CoulombScattering();
185  ss->SetEmModel(ssm);
186  ss->SetMinKinEnergy(highEnergyLimit);
187  ssm->SetLowEnergyLimit(highEnergyLimit);
188  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
189 
190  // ionisation
191  G4eIonisation* eioni = new G4eIonisation();
192 
193  // bremsstrahlung
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);
202 
203  G4ePairProduction* ee = new G4ePairProduction();
204 
205  // register processes
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);
211 
212  // e+
213  particle = G4Positron::Positron();
214 
215  // multiple scattering
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);
223 
224  // msc for HCAL using the Urban model
225  if (nullptr != aRegion) {
226  G4UrbanMscModel* msc4 = new G4UrbanMscModel();
227  msc4->SetHighEnergyLimit(highEnergyLimit);
228  msc4->SetRangeFactor(fRangeFactor);
229  msc4->SetGeomFactor(fGeomFactor);
230  msc4->SetSafetyFactor(fSafetyFactor);
231  msc4->SetLambdaLimit(fLambdaLimit);
232  msc4->SetStepLimitType(fStepLimitType);
233  msc4->SetLocked(true);
234  msc->AddEmModel(-1, msc4, aRegion);
235  }
236 
237  // msc GS with Mott-correction
238  if (nullptr != bRegion) {
239  G4GoudsmitSaundersonMscModel* msc3 = new G4GoudsmitSaundersonMscModel();
240  msc3->SetHighEnergyLimit(highEnergyLimit);
241  msc3->SetRangeFactor(0.08);
242  msc3->SetSkin(3.);
243  msc3->SetStepLimitType(fUseSafetyPlus);
244  msc3->SetLocked(true);
245  msc->AddEmModel(-1, msc3, bRegion);
246  }
247 
248  // single scattering
249  ssm = new G4eCoulombScatteringModel();
250  ss = new G4CoulombScattering();
251  ss->SetEmModel(ssm);
252  ss->SetMinKinEnergy(highEnergyLimit);
253  ssm->SetLowEnergyLimit(highEnergyLimit);
254  ssm->SetActivationLowEnergyLimit(highEnergyLimit);
255 
256  // ionisation
257  eioni = new G4eIonisation();
258 
259  // bremsstrahlung
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);
268 
269  // register processes
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);
276 
277  // generic ion
278  particle = G4GenericIon::GenericIon();
279  G4ionIonisation* ionIoni = new G4ionIonisation();
280  ph->RegisterProcess(hmsc, particle);
281  ph->RegisterProcess(ionIoni, particle);
282 
283  // muons, hadrons, ions
284  G4EmBuilder::ConstructCharged(hmsc, pnuc);
285 }

References conv, fwrapper::cs, nanoDQM_cff::Electron, fGeomFactor, fLambdaLimit, fRangeFactor, fSafetyFactor, fStepLimitType, GeV, contentValuesCheck::ss, and HLT_FULL_cff::verboseLevel.

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

HLT_FULL_cff.verboseLevel
verboseLevel
Definition: HLT_FULL_cff.py:8514
fwrapper::cs
unique_ptr< ClusterSequence > cs
Definition: fastjetfortran_madfks.cc:47
CMSEmStandardPhysicsXS::fSafetyFactor
G4double fSafetyFactor
Definition: CMSEmStandardPhysicsXS.h:30
contentValuesCheck.ss
ss
Definition: contentValuesCheck.py:33
CMSEmStandardPhysicsXS::fStepLimitType
G4MscStepLimitType fStepLimitType
Definition: CMSEmStandardPhysicsXS.h:32
AlCaHLTBitMon_ParallelJobs.p
def p
Definition: AlCaHLTBitMon_ParallelJobs.py:153
CMSEmStandardPhysicsXS::fGeomFactor
G4double fGeomFactor
Definition: CMSEmStandardPhysicsXS.h:29
GeV
const double GeV
Definition: MathUtil.h:16
AlCaHLTBitMon_QueryRunRegistry.string
string string
Definition: AlCaHLTBitMon_QueryRunRegistry.py:256
nanoDQM_cff.Electron
Electron
Definition: nanoDQM_cff.py:95
edm::LogVerbatim
Log< level::Info, true > LogVerbatim
Definition: MessageLogger.h:128
CMSEmStandardPhysicsXS::fRangeFactor
G4double fRangeFactor
Definition: CMSEmStandardPhysicsXS.h:28
CMSEmStandardPhysicsXS::fLambdaLimit
G4double fLambdaLimit
Definition: CMSEmStandardPhysicsXS.h:31
conv
EPOS::IO_EPOS conv
Definition: ReggeGribovPartonMCHadronizer.cc:42