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

#include <G4EmLivermorePhysics.hh>

+ Inheritance diagram for G4EmLivermorePhysics:

Public Member Functions

 G4EmLivermorePhysics (G4int ver=1, const G4String &name="G4EmLivermore")
 
 ~G4EmLivermorePhysics () 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 G4EmLivermorePhysics.hh.

Constructor & Destructor Documentation

◆ G4EmLivermorePhysics()

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

Definition at line 107 of file G4EmLivermorePhysics.cc.

108 : G4VPhysicsConstructor(pname), 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);
129}
@ bElectromagnetic
@ fUseSafetyPlus
void SetMinEnergy(G4double val)
void SetLowestElectronEnergy(G4double val)
void SetStepFunctionLightIons(G4double v1, G4double v2)
void SetNumberOfBinsPerDecade(G4int val)
static G4EmParameters * Instance()
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)

◆ ~G4EmLivermorePhysics()

G4EmLivermorePhysics::~G4EmLivermorePhysics ( )
override

Definition at line 133 of file G4EmLivermorePhysics.cc.

134{}

Member Function Documentation

◆ ConstructParticle()

void G4EmLivermorePhysics::ConstructParticle ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 138 of file G4EmLivermorePhysics.cc.

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

◆ ConstructProcess()

void G4EmLivermorePhysics::ConstructProcess ( )
overridevirtual

Implements G4VPhysicsConstructor.

Reimplemented in G4EmLivermorePolarizedPhysics.

Definition at line 146 of file G4EmLivermorePhysics.cc.

147{
148 if(verbose > 1) {
149 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
150 }
153
154 // processes used by several particles
156 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
157
158 // high energy limit for e+- scattering models
160 G4double livEnergyLimit = 1*GeV;
161
162 // nuclear stopping
163 G4double nielEnergyLimit = G4EmParameters::Instance()->MaxNIELEnergy();
165 pnuc->SetMaxKinEnergy(nielEnergyLimit);
166
167 // Add Livermore EM Processes
168
169 // Add gamma EM Processes
171
172 // photoelectric effect - Livermore model only
173 G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
174 thePhotoElectricEffect->SetEmModel(new G4LivermorePhotoElectricModel());
175
176 // Compton scattering - Livermore model only
177 G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
178 theComptonScattering->SetEmModel(new G4KleinNishinaModel());
179 G4VEmModel* comptLiv = new G4LivermoreComptonModel();
180 comptLiv->SetHighEnergyLimit(livEnergyLimit);
181 theComptonScattering->AddEmModel(0, comptLiv);
182
183 // gamma conversion
184 G4GammaConversion* theGammaConversion = new G4GammaConversion();
186 theGammaConversion->SetEmModel(convLiv);
187
188 // default Rayleigh scattering is Livermore
189 G4RayleighScattering* theRayleigh = new G4RayleighScattering();
190
191 ph->RegisterProcess(thePhotoElectricEffect, particle);
192 ph->RegisterProcess(theComptonScattering, particle);
193 ph->RegisterProcess(theGammaConversion, particle);
194 ph->RegisterProcess(theRayleigh, particle);
195
196 // e-
197 particle = G4Electron::Electron();
198
199 // multiple and single scattering
203 msc1->SetHighEnergyLimit(highEnergyLimit);
204 msc2->SetLowEnergyLimit(highEnergyLimit);
205 msc->SetEmModel(msc1);
206 msc->SetEmModel(msc2);
207
210 ss->SetEmModel(ssm);
211 ss->SetMinKinEnergy(highEnergyLimit);
212 ssm->SetLowEnergyLimit(highEnergyLimit);
213 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
214
215 // Ionisation - Livermore should be used only for low energies
216 G4eIonisation* eioni = new G4eIonisation();
217 G4VEmModel* theIoniLiv = new G4LivermoreIonisationModel();
218 theIoniLiv->SetHighEnergyLimit(0.1*CLHEP::MeV);
219 eioni->AddEmModel(0, theIoniLiv, new G4UniversalFluctuation() );
220
221 // Bremsstrahlung
227 brem->SetEmModel(br1);
228 brem->SetEmModel(br2);
229 br1->SetHighEnergyLimit(GeV);
230
231 // register processes
232 ph->RegisterProcess(msc, particle);
233 ph->RegisterProcess(eioni, particle);
234 ph->RegisterProcess(brem, particle);
235 ph->RegisterProcess(ee, particle);
236 ph->RegisterProcess(ss, particle);
237
238 // e+
239 particle = G4Positron::Positron();
240
241 // multiple and single scattering
242 msc = new G4eMultipleScattering;
243 msc1 = new G4GoudsmitSaundersonMscModel();
244 msc2 = new G4WentzelVIModel();
245 msc1->SetHighEnergyLimit(highEnergyLimit);
246 msc2->SetLowEnergyLimit(highEnergyLimit);
247 msc->SetEmModel(msc1);
248 msc->SetEmModel(msc2);
249
250 ssm = new G4eCoulombScatteringModel();
251 ss = new G4CoulombScattering();
252 ss->SetEmModel(ssm);
253 ss->SetMinKinEnergy(highEnergyLimit);
254 ssm->SetLowEnergyLimit(highEnergyLimit);
255 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
256
257 // ionisation
258 eioni = new G4eIonisation();
259
260 // Bremsstrahlung from standard
261 brem = new G4eBremsstrahlung();
262 br1 = new G4SeltzerBergerModel();
263 br2 = new G4eBremsstrahlungRelModel();
266 brem->SetEmModel(br1);
267 brem->SetEmModel(br2);
268 br1->SetHighEnergyLimit(GeV);
269
270 // register processes
271 ph->RegisterProcess(msc, particle);
272 ph->RegisterProcess(eioni, particle);
273 ph->RegisterProcess(brem, particle);
274 ph->RegisterProcess(ee, particle);
275 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
276 ph->RegisterProcess(ss, particle);
277
278 // generic ion
279 particle = G4GenericIon::GenericIon();
280 G4ionIonisation* ionIoni = new G4ionIonisation();
282 ph->RegisterProcess(hmsc, particle);
283 ph->RegisterProcess(ionIoni, particle);
284 ph->RegisterProcess(pnuc, particle);
285
286 // muons, hadrons, ions
288
289 // extra configuration
291
292}
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 SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:618
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

Referenced by G4EmLivermorePolarizedPhysics::ConstructProcess().


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