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

#include <G4EmLowEPPhysics.hh>

+ Inheritance diagram for G4EmLowEPPhysics:

Public Member Functions

 G4EmLowEPPhysics (G4int ver=1)
 
 G4EmLowEPPhysics (G4int ver, const G4String &name)
 
virtual ~G4EmLowEPPhysics ()
 
virtual void ConstructParticle ()
 
virtual void ConstructProcess ()
 
- 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
 

Additional Inherited Members

- Protected Member Functions inherited from G4VPhysicsConstructor
G4bool RegisterProcess (G4VProcess *process, G4ParticleDefinition *particle)
 
- Protected Attributes inherited from G4VPhysicsConstructor
G4int verboseLevel
 
G4String namePhysics
 
G4int typePhysics
 
G4ParticleTabletheParticleTable
 
G4ParticleTable::G4PTblDicIteratortheParticleIterator
 
G4PhysicsListHelperthePLHelper
 

Detailed Description

Definition at line 35 of file G4EmLowEPPhysics.hh.

Constructor & Destructor Documentation

◆ G4EmLowEPPhysics() [1/2]

G4EmLowEPPhysics::G4EmLowEPPhysics ( G4int  ver = 1)

Definition at line 126 of file G4EmLowEPPhysics.cc.

127 : G4VPhysicsConstructor("G4EmLowEPPhysics"), verbose(ver)
128{
131}
@ bElectromagnetic
static G4LossTableManager * Instance()

◆ G4EmLowEPPhysics() [2/2]

G4EmLowEPPhysics::G4EmLowEPPhysics ( G4int  ver,
const G4String name 
)

Definition at line 135 of file G4EmLowEPPhysics.cc.

136 : G4VPhysicsConstructor("G4EmLowEPPhysics"), verbose(ver)
137{
140}

◆ ~G4EmLowEPPhysics()

G4EmLowEPPhysics::~G4EmLowEPPhysics ( )
virtual

Definition at line 144 of file G4EmLowEPPhysics.cc.

145{}

Member Function Documentation

◆ ConstructParticle()

void G4EmLowEPPhysics::ConstructParticle ( )
virtual

Implements G4VPhysicsConstructor.

Definition at line 149 of file G4EmLowEPPhysics.cc.

150{
151// gamma
153
154// leptons
159
160// mesons
165
166// baryons
169
170// ions
173 G4He3::He3();
176}
static G4Alpha * Alpha()
Definition: G4Alpha.cc:89
static G4AntiProton * AntiProton()
Definition: G4AntiProton.cc:93
static G4Deuteron * Deuteron()
Definition: G4Deuteron.cc:94
static G4Electron * Electron()
Definition: G4Electron.cc:94
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
static G4GenericIon * GenericIonDefinition()
Definition: G4GenericIon.cc:87
static G4He3 * He3()
Definition: G4He3.cc:94
static G4KaonMinus * KaonMinusDefinition()
Definition: G4KaonMinus.cc:108
static G4KaonPlus * KaonPlusDefinition()
Definition: G4KaonPlus.cc:108
static G4MuonMinus * MuonMinus()
Definition: G4MuonMinus.cc:100
static G4MuonPlus * MuonPlus()
Definition: G4MuonPlus.cc:99
static G4PionMinus * PionMinusDefinition()
Definition: G4PionMinus.cc:93
static G4PionPlus * PionPlusDefinition()
Definition: G4PionPlus.cc:93
static G4Positron * Positron()
Definition: G4Positron.cc:94
static G4Proton * Proton()
Definition: G4Proton.cc:93
static G4Triton * Triton()
Definition: G4Triton.cc:95

◆ ConstructProcess()

void G4EmLowEPPhysics::ConstructProcess ( )
virtual

Implements G4VPhysicsConstructor.

Definition at line 180 of file G4EmLowEPPhysics.cc.

181{
183
184 // muon & hadron bremsstrahlung and pair production
193
194 // muon & hadron multiple scattering
196 mumsc->AddEmModel(0, new G4WentzelVIModel());
201
202 // Add Livermore EM Processes
204
205 while( (*theParticleIterator)() ){
206
208 G4String particleName = particle->GetParticleName();
209
210 if(verbose > 1)
211 G4cout << "### " << GetPhysicsName() << " instantiates for "
212 << particleName << G4endl;
213
214 //Applicability range for Livermore models
215 //for higher energies, the Standard models are used
216 G4double LivermoreHighEnergyLimit = GeV;
217
218 if (particleName == "gamma") {
219
220 G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
221 G4LivermorePhotoElectricModel* theLivermorePhotoElectricModel =
223 theLivermorePhotoElectricModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
224 thePhotoElectricEffect->AddEmModel(0, theLivermorePhotoElectricModel);
225 ph->RegisterProcess(thePhotoElectricEffect, particle);
226
227 G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
228 G4LowEPComptonModel* theLowEPComptonModel =
230 theLowEPComptonModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
231 theComptonScattering->AddEmModel(0, theLowEPComptonModel);
232 ph->RegisterProcess(theComptonScattering, particle);
233
234 G4GammaConversion* theGammaConversion = new G4GammaConversion();
235 G4LivermoreGammaConversionModel* theLivermoreGammaConversionModel =
237 theLivermoreGammaConversionModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
238 theGammaConversion->AddEmModel(0, theLivermoreGammaConversionModel);
239 ph->RegisterProcess(theGammaConversion, particle);
240
241 G4RayleighScattering* theRayleigh = new G4RayleighScattering();
242 G4LivermoreRayleighModel* theRayleighModel = new G4LivermoreRayleighModel();
243 theRayleighModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
244 theRayleigh->AddEmModel(0, theRayleighModel);
245 ph->RegisterProcess(theRayleigh, particle);
246
247 } else if (particleName == "e-") {
248
250 msc->AddEmModel(0, new G4UrbanMscModel95());
251 //msc->AddEmModel(0, new G4GoudsmitSaundersonMscModel());
253 ph->RegisterProcess(msc, particle);
254
255 // Ionisation
256 G4eIonisation* eIoni = new G4eIonisation();
257 G4LivermoreIonisationModel* theIoniLivermore = new
259 theIoniLivermore->SetHighEnergyLimit(0.1*MeV);
260 eIoni->AddEmModel(0, theIoniLivermore, new G4UniversalFluctuation() );
261 eIoni->SetStepFunction(0.2, 100*um); //
262 ph->RegisterProcess(eIoni, particle);
263
264 // Bremsstrahlung
266 G4LivermoreBremsstrahlungModel* theBremLivermore = new
268 theBremLivermore->SetHighEnergyLimit(25*MeV);
269 theBremLivermore->SetAngularDistribution(new G4Generator2BS());
270 eBrem->AddEmModel(0, theBremLivermore);
271 ph->RegisterProcess(eBrem, particle);
272
273 } else if (particleName == "e+") {
274
275 // Identical to G4EmStandardPhysics_option3
276
278 msc->AddEmModel(0, new G4UrbanMscModel95());
279 //msc->AddEmModel(0, new G4GoudsmitSaundersonMscModel());
281 G4eIonisation* eIoni = new G4eIonisation();
282 eIoni->SetStepFunction(0.2, 100*um);
283
284 ph->RegisterProcess(msc, particle);
285 ph->RegisterProcess(eIoni, particle);
286 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
287 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
288
289 } else if (particleName == "mu+" ||
290 particleName == "mu-" ) {
291
292 G4MuIonisation* muIoni = new G4MuIonisation();
293 muIoni->SetStepFunction(0.2, 50*um);
294
295 ph->RegisterProcess(mumsc, particle);
296 ph->RegisterProcess(muIoni, particle);
297 ph->RegisterProcess(mub, particle);
298 ph->RegisterProcess(mup, particle);
299 ph->RegisterProcess(new G4CoulombScattering(), particle);
300
301 } else if (particleName == "alpha" ||
302 particleName == "He3" ) {
303
304 // Identical to G4EmStandardPhysics_option3
305
306 G4ionIonisation* ionIoni = new G4ionIonisation();
307 ionIoni->SetStepFunction(0.1, 10*um);
308
309 ph->RegisterProcess(hmsc, particle);
310 ph->RegisterProcess(ionIoni, particle);
311 ph->RegisterProcess(new G4NuclearStopping(), particle);
312
313 } else if (particleName == "GenericIon") {
314
315 // Identical to G4EmStandardPhysics_option3
316
317 G4ionIonisation* ionIoni = new G4ionIonisation();
319 ionIoni->SetStepFunction(0.1, 1*um);
320
321 ph->RegisterProcess(hmsc, particle);
322 ph->RegisterProcess(ionIoni, particle);
323 ph->RegisterProcess(new G4NuclearStopping(), particle);
324
325 } else if (particleName == "pi+" ||
326 particleName == "pi-" ) {
327
328 G4hIonisation* hIoni = new G4hIonisation();
329 hIoni->SetStepFunction(0.2, 50*um);
330
331 ph->RegisterProcess(pimsc, particle);
332 ph->RegisterProcess(hIoni, particle);
333 ph->RegisterProcess(pib, particle);
334 ph->RegisterProcess(pip, particle);
335
336 } else if (particleName == "kaon+" ||
337 particleName == "kaon-" ) {
338
339 G4hIonisation* hIoni = new G4hIonisation();
340 hIoni->SetStepFunction(0.2, 50*um);
341
342 ph->RegisterProcess(kmsc, particle);
343 ph->RegisterProcess(hIoni, particle);
344 ph->RegisterProcess(kb, particle);
345 ph->RegisterProcess(kp, particle);
346
347 } else if (particleName == "proton" ||
348 particleName == "anti_proton") {
349
350 G4hIonisation* hIoni = new G4hIonisation();
351 hIoni->SetStepFunction(0.2, 50*um);
352
353 ph->RegisterProcess(pmsc, particle);
354 ph->RegisterProcess(hIoni, particle);
355 ph->RegisterProcess(pb, particle);
356 ph->RegisterProcess(pp, particle);
357
358 } else if (particleName == "B+" ||
359 particleName == "B-" ||
360 particleName == "D+" ||
361 particleName == "D-" ||
362 particleName == "Ds+" ||
363 particleName == "Ds-" ||
364 particleName == "anti_He3" ||
365 particleName == "anti_alpha" ||
366 particleName == "anti_deuteron" ||
367 particleName == "anti_lambda_c+" ||
368 particleName == "anti_omega-" ||
369 particleName == "anti_sigma_c+" ||
370 particleName == "anti_sigma_c++" ||
371 particleName == "anti_sigma+" ||
372 particleName == "anti_sigma-" ||
373 particleName == "anti_triton" ||
374 particleName == "anti_xi_c+" ||
375 particleName == "anti_xi-" ||
376 particleName == "deuteron" ||
377 particleName == "lambda_c+" ||
378 particleName == "omega-" ||
379 particleName == "sigma_c+" ||
380 particleName == "sigma_c++" ||
381 particleName == "sigma+" ||
382 particleName == "sigma-" ||
383 particleName == "tau+" ||
384 particleName == "tau-" ||
385 particleName == "triton" ||
386 particleName == "xi_c+" ||
387 particleName == "xi-" ) {
388
389 // Identical to G4EmStandardPhysics_option3
390
391 ph->RegisterProcess(hmsc, particle);
392 ph->RegisterProcess(new G4hIonisation(), particle);
393
394 }
395 }
396
397 // Em options
398 //
400 opt.SetVerbose(verbose);
401
402 // Multiple Coulomb scattering
403 //
404 opt.SetPolarAngleLimit(CLHEP::pi);
405
406 // Physics tables
407 //
408
409 opt.SetMinEnergy(100*eV);
410 opt.SetMaxEnergy(10*TeV);
411 opt.SetDEDXBinning(220);
412 opt.SetLambdaBinning(220);
413
414 // Ionization
415 //
416 //opt.SetSubCutoff(true);
417
418 // Deexcitation
419 //
422 de->SetFluo(true);
423}
@ fUseDistanceToBoundary
double G4double
Definition: G4Types.hh:64
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
void SetMaxEnergy(G4double val)
void SetDEDXBinning(G4int val)
void SetLambdaBinning(G4int val)
void SetPolarAngleLimit(G4double val)
void SetVerbose(G4int val, const G4String &name="all")
void SetMinEnergy(G4double val)
void SetAtomDeexcitation(G4VAtomDeexcitation *)
const G4String & GetParticleName() const
G4bool RegisterProcess(G4VProcess *process, G4ParticleDefinition *particle)
static G4PhysicsListHelper * GetPhysicsListHelper()
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:515
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=0)
void SetEmModel(G4VEmModel *, G4int index=1)
void AddEmModel(G4int, G4VEmModel *, G4VEmFluctuationModel *fluc=0, const G4Region *region=0)
void SetStepFunction(G4double v1, G4double v2)
void AddEmModel(G4int order, G4VEmModel *, const G4Region *region=0)
void SetStepLimitType(G4MscStepLimitType val)
const G4String & GetPhysicsName() const
G4ParticleTable::G4PTblDicIterator * theParticleIterator

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