Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4EmPenelopePhysics Class Reference

#include <G4EmPenelopePhysics.hh>

+ Inheritance diagram for G4EmPenelopePhysics:

Public Member Functions

 G4EmPenelopePhysics (G4int ver=1, const G4String &name="")
 
 ~G4EmPenelopePhysics () override
 
void ConstructParticle () override
 
void ConstructProcess () override
 
- Public Member Functions inherited from G4VPhysicsConstructor
 G4VPhysicsConstructor (const G4String &="")
 
 G4VPhysicsConstructor (const G4String &name, G4int physics_type)
 
virtual ~G4VPhysicsConstructor ()
 
virtual void ConstructParticle ()=0
 
virtual void ConstructProcess ()=0
 
void SetPhysicsName (const G4String &="")
 
const G4StringGetPhysicsName () const
 
void SetPhysicsType (G4int)
 
G4int GetPhysicsType () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
G4int GetInstanceID () const
 
virtual void TerminateWorker ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VPhysicsConstructor
static const G4VPCManagerGetSubInstanceManager ()
 
- Protected Types inherited from G4VPhysicsConstructor
using PhysicsBuilder_V = G4VPCData::PhysicsBuilders_V
 
- Protected Member Functions inherited from G4VPhysicsConstructor
G4bool RegisterProcess (G4VProcess *process, G4ParticleDefinition *particle)
 
G4ParticleTable::G4PTblDicIteratorGetParticleIterator () const
 
PhysicsBuilder_V GetBuilders () const
 
void AddBuilder (G4PhysicsBuilderInterface *bld)
 
- Protected Attributes inherited from G4VPhysicsConstructor
G4int verboseLevel
 
G4String namePhysics
 
G4int typePhysics
 
G4ParticleTabletheParticleTable
 
G4int g4vpcInstanceID
 
- Static Protected Attributes inherited from G4VPhysicsConstructor
static G4RUN_DLL G4VPCManager subInstanceManager
 

Detailed Description

Definition at line 35 of file G4EmPenelopePhysics.hh.

Constructor & Destructor Documentation

◆ G4EmPenelopePhysics()

G4EmPenelopePhysics::G4EmPenelopePhysics ( G4int  ver = 1,
const G4String name = "" 
)
explicit

Definition at line 107 of file G4EmPenelopePhysics.cc.

108 : G4VPhysicsConstructor("G4EmPenelope"), verbose(ver)
109{
111 param->SetDefaults();
112 param->SetVerbose(verbose);
113 param->SetMinEnergy(100*CLHEP::eV);
114 param->SetLowestElectronEnergy(100*CLHEP::eV);
115 param->SetNumberOfBinsPerDecade(20);
117 param->SetStepFunction(0.2, 10*CLHEP::um);
118 param->SetStepFunctionMuHad(0.1, 50*CLHEP::um);
119 param->SetStepFunctionLightIons(0.1, 20*CLHEP::um);
120 param->SetStepFunctionIons(0.1, 1*CLHEP::um);
121 param->SetUseMottCorrection(true);
123 param->SetMscSkin(3);
124 param->SetMscRangeFactor(0.08);
125 param->SetMuHadLateralDisplacement(true);
126 param->SetFluo(true);
127 param->SetMaxNIELEnergy(1*CLHEP::MeV);
128 param->SetPIXEElectronCrossSectionModel("Penelope");
130}
@ bElectromagnetic
@ fUseSafetyPlus
void SetMinEnergy(G4double val)
void SetLowestElectronEnergy(G4double val)
void SetStepFunctionLightIons(G4double v1, G4double v2)
void SetNumberOfBinsPerDecade(G4int val)
static G4EmParameters * Instance()
void SetPIXEElectronCrossSectionModel(const G4String &)
void SetMuHadLateralDisplacement(G4bool val)
void ActivateAngularGeneratorForIonisation(G4bool val)
void SetStepFunction(G4double v1, G4double v2)
void SetFluo(G4bool val)
void SetStepFunctionMuHad(G4double v1, G4double v2)
void SetVerbose(G4int val)
void SetMscSkin(G4double val)
void SetMaxNIELEnergy(G4double val)
void SetStepFunctionIons(G4double v1, G4double v2)
void SetMscStepLimitType(G4MscStepLimitType val)
void SetUseMottCorrection(G4bool val)
void SetMscRangeFactor(G4double val)

◆ ~G4EmPenelopePhysics()

G4EmPenelopePhysics::~G4EmPenelopePhysics ( )
override

Definition at line 134 of file G4EmPenelopePhysics.cc.

135{}

Member Function Documentation

◆ ConstructParticle()

void G4EmPenelopePhysics::ConstructParticle ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 139 of file G4EmPenelopePhysics.cc.

140{
141 // minimal set of particles for EM physics
143}
static void ConstructMinimalEmSet()
Definition: G4EmBuilder.cc:234

◆ ConstructProcess()

void G4EmPenelopePhysics::ConstructProcess ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 147 of file G4EmPenelopePhysics.cc.

148{
149 if(verbose > 1) {
150 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
151 }
154
155 // processes used by several particles
157 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
158
159 // high energy limit for e+- scattering models
160 G4double highEnergyLimit = G4EmParameters::Instance()->MscEnergyLimit();
161
162 // nuclear stopping
163 G4double nielEnergyLimit = G4EmParameters::Instance()->MaxNIELEnergy();
165 pnuc->SetMaxKinEnergy(nielEnergyLimit);
166
167 //Applicability range for Penelope models
168 //for higher energies, the Standard models are used
169 G4double PenelopeHighEnergyLimit = 1.0*CLHEP::GeV;
170
171 // Add gamma EM Processes
173
174 //Photo-electric effect
175 G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
176 thePhotoElectricEffect->SetEmModel(new G4PEEffectFluoModel());
178 thePEPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
179 thePhotoElectricEffect->AddEmModel(0, thePEPenelopeModel);
180
181 //Compton scattering
182 G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
183 G4PenelopeComptonModel* theComptonPenelopeModel = new G4PenelopeComptonModel();
184 theComptonScattering->SetEmModel(new G4KleinNishinaModel());
185 theComptonPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
186 theComptonScattering->AddEmModel(0, theComptonPenelopeModel);
187
188 //Gamma conversion
189 G4GammaConversion* theGammaConversion = new G4GammaConversion();
191 theGCPenelopeModel->SetHighEnergyLimit(PenelopeHighEnergyLimit);
192 theGammaConversion->AddEmModel(0, theGCPenelopeModel);
193
194 //Rayleigh scattering
195 G4RayleighScattering* theRayleigh = new G4RayleighScattering();
196 G4PenelopeRayleighModel* theRayleighPenelopeModel = new G4PenelopeRayleighModel();
197 theRayleigh->SetEmModel(theRayleighPenelopeModel);
198
199 ph->RegisterProcess(thePhotoElectricEffect, particle);
200 ph->RegisterProcess(theComptonScattering, particle);
201 ph->RegisterProcess(theGammaConversion, particle);
202 ph->RegisterProcess(theRayleigh, particle);
203
204 // e-
205 particle = G4Electron::Electron();
206
207 // multiple scattering
211 msc1->SetHighEnergyLimit(highEnergyLimit);
212 msc2->SetLowEnergyLimit(highEnergyLimit);
213 msc->SetEmModel(msc1);
214 msc->SetEmModel(msc2);
215
218 ss->SetEmModel(ssm);
219 ss->SetMinKinEnergy(highEnergyLimit);
220 ssm->SetLowEnergyLimit(highEnergyLimit);
221 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
222
223 //Ionisation
224 G4eIonisation* eioni = new G4eIonisation();
226 theIoniPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
227 eioni->AddEmModel(0, theIoniPenelope, new G4UniversalFluctuation());
228
229 //Bremsstrahlung
232 theBremPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
233 brem->SetEmModel(theBremPenelope);
234
235 // register processes
236 ph->RegisterProcess(msc, particle);
237 ph->RegisterProcess(eioni, particle);
238 ph->RegisterProcess(brem, particle);
239 ph->RegisterProcess(ee, particle);
240 ph->RegisterProcess(ss, particle);
241
242 // e+
243 particle = G4Positron::Positron();
244
245 // multiple scattering
246 msc = new G4eMultipleScattering;
247 msc1 = new G4GoudsmitSaundersonMscModel();
248 msc2 = new G4WentzelVIModel();
249 msc1->SetHighEnergyLimit(highEnergyLimit);
250 msc2->SetLowEnergyLimit(highEnergyLimit);
251 msc->SetEmModel(msc1);
252 msc->SetEmModel(msc2);
253
254 ssm = new G4eCoulombScatteringModel();
255 ss = new G4CoulombScattering();
256 ss->SetEmModel(ssm);
257 ss->SetMinKinEnergy(highEnergyLimit);
258 ssm->SetLowEnergyLimit(highEnergyLimit);
259 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
260
261 //Ionisation
262 eioni = new G4eIonisation();
263 theIoniPenelope = new G4PenelopeIonisationModel();
264 theIoniPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
265 eioni->AddEmModel(0,theIoniPenelope, new G4UniversalFluctuation());
266
267 //Bremsstrahlung
268 brem = new G4eBremsstrahlung();
269 theBremPenelope = new G4PenelopeBremsstrahlungModel();
270 theBremPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
271 brem->SetEmModel(theBremPenelope);
272
273 //Annihilation
276 theAnnPenelope->SetHighEnergyLimit(PenelopeHighEnergyLimit);
277 anni->AddEmModel(0, theAnnPenelope);
278
279 // register processes
280 ph->RegisterProcess(msc, particle);
281 ph->RegisterProcess(eioni, particle);
282 ph->RegisterProcess(brem, particle);
283 ph->RegisterProcess(ee, particle);
284 ph->RegisterProcess(anni, particle);
285 ph->RegisterProcess(ss, particle);
286
287 // generic ion
288 particle = G4GenericIon::GenericIon();
289 G4ionIonisation* ionIoni = new G4ionIonisation();
291 ph->RegisterProcess(hmsc, particle);
292 ph->RegisterProcess(ionIoni, particle);
293 ph->RegisterProcess(pnuc, particle);
294
295 // muons, hadrons, ions
297
298 // extra configuration
300}
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
static G4Electron * Electron()
Definition: G4Electron.cc:93
static void ConstructCharged(G4hMultipleScattering *hmsc, G4NuclearStopping *nucStopping, G4bool isWVI=true)
Definition: G4EmBuilder.cc:170
static void PrepareEMPhysics()
Definition: G4EmBuilder.cc:259
G4double MaxNIELEnergy() const
G4double MscEnergyLimit() const
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
static G4GenericIon * GenericIon()
Definition: G4GenericIon.cc:92
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4PhysicsListHelper * GetPhysicsListHelper()
static G4Positron * Positron()
Definition: G4Positron.cc:93
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757
void SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:778
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:764
void SetMinKinEnergy(G4double e)
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetMaxKinEnergy(G4double e)
void SetEmModel(G4VEmModel *, G4int index=0)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=nullptr)
void SetEmModel(G4VMscModel *, size_t index=0)
const G4String & GetPhysicsName() const

The documentation for this class was generated from the following files: