Geant4 11.3.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 ()
 
void SetPhysicsName (const G4String &="")
 
const G4StringGetPhysicsName () const
 
void SetPhysicsType (G4int)
 
G4int GetPhysicsType () const
 
G4int GetInstanceID () const
 
virtual void TerminateWorker ()
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 

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 = 0
 
G4String namePhysics = ""
 
G4int typePhysics = 0
 
G4ParticleTabletheParticleTable = nullptr
 
G4int g4vpcInstanceID = 0
 
- 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 112 of file G4EmLivermorePhysics.cc.

113 : G4VPhysicsConstructor(pname)
114{
115 SetVerboseLevel(ver);
116 G4EmParameters* param = G4EmParameters::Instance();
117 param->SetDefaults();
118 param->SetVerbose(ver);
119 param->SetMinEnergy(100*CLHEP::eV);
120 param->SetLowestElectronEnergy(100*CLHEP::eV);
121 param->SetNumberOfBinsPerDecade(20);
123 param->SetStepFunction(0.2, 10*CLHEP::um);
124 param->SetStepFunctionMuHad(0.1, 50*CLHEP::um);
125 param->SetStepFunctionLightIons(0.1, 20*CLHEP::um);
126 param->SetStepFunctionIons(0.1, 1*CLHEP::um);
127 param->SetUseMottCorrection(true);
129 param->SetMscSkin(3);
130 param->SetMscRangeFactor(0.08);
131 param->SetMuHadLateralDisplacement(true);
132 param->SetFluo(true);
133 param->SetUseICRU90Data(true);
135 param->SetMaxNIELEnergy(1*CLHEP::MeV);
138}
@ bElectromagnetic
@ fUrbanFluctuation
@ fAllisonPositronium
@ fUseSafetyPlus
void SetMinEnergy(G4double val)
void SetLowestElectronEnergy(G4double val)
void SetStepFunctionLightIons(G4double v1, G4double v2)
void SetNumberOfBinsPerDecade(G4int val)
static G4EmParameters * Instance()
void SetPositronAtRestModelType(G4PositronAtRestModelType val)
void SetMuHadLateralDisplacement(G4bool val)
void ActivateAngularGeneratorForIonisation(G4bool val)
void SetStepFunction(G4double v1, G4double v2)
void SetFluo(G4bool val)
void SetFluctuationType(G4EmFluctuationType 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 SetUseICRU90Data(G4bool val)
void SetUseMottCorrection(G4bool val)
void SetMscRangeFactor(G4double val)
G4VPhysicsConstructor(const G4String &="")
void SetVerboseLevel(G4int value)

Referenced by G4EmLivermorePolarizedPhysics::G4EmLivermorePolarizedPhysics().

◆ ~G4EmLivermorePhysics()

G4EmLivermorePhysics::~G4EmLivermorePhysics ( )
overridedefault

Member Function Documentation

◆ ConstructParticle()

void G4EmLivermorePhysics::ConstructParticle ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 146 of file G4EmLivermorePhysics.cc.

147{
148 // minimal set of particles for EM physics
150}
static void ConstructMinimalEmSet()

◆ ConstructProcess()

void G4EmLivermorePhysics::ConstructProcess ( )
overridevirtual

Implements G4VPhysicsConstructor.

Definition at line 154 of file G4EmLivermorePhysics.cc.

155{
156 if(verboseLevel > 1) {
157 G4cout << "### " << GetPhysicsName() << " Construct Processes " << G4endl;
158 }
160 G4PhysicsListHelper* ph = G4PhysicsListHelper::GetPhysicsListHelper();
161 G4EmParameters* param = G4EmParameters::Instance();
162
163 // processes used by several particles
164 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
165
166 // high energy limit for e+- scattering models
167 G4double highEnergyLimit= param->MscEnergyLimit();
168 G4double livEnergyLimit = 1*CLHEP::GeV;
169
170 // nuclear stopping
171 G4double nielEnergyLimit = param->MaxNIELEnergy();
172 G4NuclearStopping* pnuc = nullptr;
173 if(nielEnergyLimit > 0.0) {
174 pnuc = new G4NuclearStopping();
175 pnuc->SetMaxKinEnergy(nielEnergyLimit);
176 }
177
178 // Add Livermore EM Processes
179
180 // Add gamma EM Processes
181 G4ParticleDefinition* particle = G4Gamma::Gamma();
182 G4bool polar = param->EnablePolarisation();
183
184 // photoelectric effect - Livermore model only
185 G4PhotoElectricEffect* pe = new G4PhotoElectricEffect();
186 G4VEmModel* peModel = new G4LivermorePhotoElectricModel();
187 pe->SetEmModel(peModel);
188 if(polar) {
189 peModel->SetAngularDistribution(new G4PhotoElectricAngularGeneratorPolarized());
190 }
191
192 // Compton scattering - Livermore model only
193 G4ComptonScattering* cs = new G4ComptonScattering;
194 cs->SetEmModel(new G4KleinNishinaModel());
195 G4VEmModel* cModel = nullptr;
196 if(polar) {
197 cModel = new G4LivermorePolarizedComptonModel();
198 cModel->SetHighEnergyLimit(20*CLHEP::MeV);
199 } else {
200 cModel = new G4LivermoreComptonModel();
201 cModel->SetHighEnergyLimit(livEnergyLimit);
202 }
203 cs->AddEmModel(0, cModel);
204
205 // gamma conversion
206 G4GammaConversion* gc = new G4GammaConversion();
207 G4VEmModel* convLiv = new G4LivermoreGammaConversion5DModel();
208 gc->SetEmModel(convLiv);
209
210 // Rayleigh scattering
211 G4RayleighScattering* rl = new G4RayleighScattering();
212 if(polar) {
213 rl->SetEmModel(new G4LivermorePolarizedRayleighModel());
214 }
215
216 ph->RegisterProcess(pe, particle);
217 ph->RegisterProcess(cs, particle);
218 ph->RegisterProcess(gc, particle);
219 ph->RegisterProcess(rl, particle);
220
221 // e-
222 particle = G4Electron::Electron();
223
224 // multiple and single scattering
225 G4GoudsmitSaundersonMscModel* msc1 = new G4GoudsmitSaundersonMscModel();
226 G4WentzelVIModel* msc2 = new G4WentzelVIModel();
227 msc1->SetHighEnergyLimit(highEnergyLimit);
228 msc2->SetLowEnergyLimit(highEnergyLimit);
229 G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle);
230
231 G4eCoulombScatteringModel* ssm = new G4eCoulombScatteringModel();
232 G4CoulombScattering* ss = new G4CoulombScattering();
233 ss->SetEmModel(ssm);
234 ss->SetMinKinEnergy(highEnergyLimit);
235 ssm->SetLowEnergyLimit(highEnergyLimit);
236 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
237
238 // Ionisation - Livermore should be used only for low energies
239 G4eIonisation* eioni = new G4eIonisation();
241 G4VEmModel* theIoniLiv = new G4LivermoreIonisationModel();
242 theIoniLiv->SetHighEnergyLimit(0.1*CLHEP::MeV);
243 eioni->AddEmModel(0, theIoniLiv);
244
245 // Bremsstrahlung
246 G4eBremsstrahlung* brem = new G4eBremsstrahlung();
247 G4SeltzerBergerModel* br1 = new G4SeltzerBergerModel();
248 G4eBremsstrahlungRelModel* br2 = new G4eBremsstrahlungRelModel();
249 br1->SetAngularDistribution(new G4Generator2BS());
250 br2->SetAngularDistribution(new G4Generator2BS());
251 brem->SetEmModel(br1);
252 brem->SetEmModel(br2);
253 br1->SetHighEnergyLimit(GeV);
254
255 G4ePairProduction* ee = new G4ePairProduction();
256
257 // register processes
258 ph->RegisterProcess(eioni, particle);
259 ph->RegisterProcess(brem, particle);
260 ph->RegisterProcess(ee, particle);
261 ph->RegisterProcess(ss, particle);
262
263 // e+
264 particle = G4Positron::Positron();
265
266 // multiple and single scattering
267 msc1 = new G4GoudsmitSaundersonMscModel();
268 msc2 = new G4WentzelVIModel();
269 msc1->SetHighEnergyLimit(highEnergyLimit);
270 msc2->SetLowEnergyLimit(highEnergyLimit);
271 G4EmBuilder::ConstructElectronMscProcess(msc1, msc2, particle);
272
273 ssm = new G4eCoulombScatteringModel();
274 ss = new G4CoulombScattering();
275 ss->SetEmModel(ssm);
276 ss->SetMinKinEnergy(highEnergyLimit);
277 ssm->SetLowEnergyLimit(highEnergyLimit);
278 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
279
280 // ionisation
281 eioni = new G4eIonisation();
282
283 // Bremsstrahlung from standard
284 brem = new G4eBremsstrahlung();
285 br1 = new G4SeltzerBergerModel();
286 br2 = new G4eBremsstrahlungRelModel();
287 br1->SetAngularDistribution(new G4Generator2BS());
288 br2->SetAngularDistribution(new G4Generator2BS());
289 brem->SetEmModel(br1);
290 brem->SetEmModel(br2);
291 br1->SetHighEnergyLimit(GeV);
292
293 // annihilation
294 auto anni = new G4eplusAnnihilation();
295 if (param->Use3GammaAnnihilationOnFly()) {
296 anni->SetEmModel(new G4eplusTo2or3GammaModel());
297 }
298
299 // register processes
300 ph->RegisterProcess(eioni, particle);
301 ph->RegisterProcess(brem, particle);
302 ph->RegisterProcess(ee, particle);
303 ph->RegisterProcess(anni, particle);
304 ph->RegisterProcess(ss, particle);
305
306 // generic ion
307 particle = G4GenericIon::GenericIon();
308 G4ionIonisation* ionIoni = new G4ionIonisation();
309 ionIoni->SetEmModel(new G4LindhardSorensenIonModel());
310 ph->RegisterProcess(hmsc, particle);
311 ph->RegisterProcess(ionIoni, particle);
312 if(nullptr != pnuc) { ph->RegisterProcess(pnuc, particle); }
313
314 // muons, hadrons, ions
316
317 // extra configuration
318 G4EmModelActivator mact(GetPhysicsName());
319
320}
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
static G4Electron * Electron()
Definition G4Electron.cc:91
static void ConstructCharged(G4hMultipleScattering *hmsc, G4NuclearStopping *nucStopping, G4bool isWVI=true)
static void ConstructElectronMscProcess(G4VMscModel *msc1, G4VMscModel *msc2, G4ParticleDefinition *particle)
static void PrepareEMPhysics()
G4bool EnablePolarisation() const
G4double MaxNIELEnergy() const
G4double MscEnergyLimit() const
G4bool Use3GammaAnnihilationOnFly() const
static G4VEmFluctuationModel * ModelOfFluctuations(G4bool isIon=false)
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
static G4GenericIon * GenericIon()
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4PhysicsListHelper * GetPhysicsListHelper()
static G4Positron * Positron()
Definition G4Positron.cc:90
void SetHighEnergyLimit(G4double)
void SetActivationLowEnergyLimit(G4double)
void SetLowEnergyLimit(G4double)
void SetAngularDistribution(G4VEmAngularDistribution *)
void SetMinKinEnergy(G4double e)
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
void SetEmModel(G4VEmModel *, G4int index=0)
void SetMaxKinEnergy(G4double e)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=nullptr, const G4Region *region=nullptr)
void SetFluctModel(G4VEmFluctuationModel *)
void SetEmModel(G4VEmModel *, G4int index=0)
const G4String & GetPhysicsName() const

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