Geant4 9.6.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)
 
 G4EmLivermorePhysics (G4int ver, const G4String &name)
 
virtual ~G4EmLivermorePhysics ()
 
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 36 of file G4EmLivermorePhysics.hh.

Constructor & Destructor Documentation

◆ G4EmLivermorePhysics() [1/2]

G4EmLivermorePhysics::G4EmLivermorePhysics ( G4int  ver = 1)

Definition at line 128 of file G4EmLivermorePhysics.cc.

129 : G4VPhysicsConstructor("G4EmLivermorePhysics"), verbose(ver)
130{
133}
@ bElectromagnetic
static G4LossTableManager * Instance()

◆ G4EmLivermorePhysics() [2/2]

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

Definition at line 137 of file G4EmLivermorePhysics.cc.

138 : G4VPhysicsConstructor("G4EmLivermorePhysics"), verbose(ver)
139{
142}

◆ ~G4EmLivermorePhysics()

G4EmLivermorePhysics::~G4EmLivermorePhysics ( )
virtual

Definition at line 146 of file G4EmLivermorePhysics.cc.

147{}

Member Function Documentation

◆ ConstructParticle()

void G4EmLivermorePhysics::ConstructParticle ( )
virtual

Implements G4VPhysicsConstructor.

Definition at line 151 of file G4EmLivermorePhysics.cc.

152{
153// gamma
155
156// leptons
161
162// mesons
167
168// baryons
171
172// ions
175 G4He3::He3();
178}
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 G4EmLivermorePhysics::ConstructProcess ( )
virtual

Implements G4VPhysicsConstructor.

Definition at line 182 of file G4EmLivermorePhysics.cc.

183{
185
186 // muon & hadron bremsstrahlung and pair production
195
196 // muon & hadron multiple scattering
198 mumsc->AddEmModel(0, new G4WentzelVIModel());
200 pimsc->AddEmModel(0, new G4WentzelVIModel());
202 kmsc->AddEmModel(0, new G4WentzelVIModel());
204 pmsc->AddEmModel(0, new G4WentzelVIModel());
205 G4hMultipleScattering* hmsc = new G4hMultipleScattering("ionmsc");
206
207 // high energy limit for e+- scattering models
208 G4double highEnergyLimit = 100*MeV;
209
210 // nuclear stopping
211 G4NuclearStopping* ionnuc = new G4NuclearStopping();
213
214 // Add Livermore EM Processes
216
217 while( (*theParticleIterator)() ){
218
220 G4String particleName = particle->GetParticleName();
221
222 if(verbose > 1)
223 G4cout << "### " << GetPhysicsName() << " instantiates for "
224 << particleName << G4endl;
225
226 //Applicability range for Livermore models
227 //for higher energies, the Standard models are used
228 G4double LivermoreHighEnergyLimit = GeV;
229
230 if (particleName == "gamma") {
231
232 // Photoelectric effect - define low-energy model
233 G4PhotoElectricEffect* thePhotoElectricEffect = new G4PhotoElectricEffect();
234 G4LivermorePhotoElectricModel* theLivermorePhotoElectricModel =
236 thePhotoElectricEffect->SetEmModel(theLivermorePhotoElectricModel);
237 ph->RegisterProcess(thePhotoElectricEffect, particle);
238
239 // Compton scattering - define low-energy model
240 G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
241 G4LivermoreComptonModel* theLivermoreComptonModel =
243 theLivermoreComptonModel->SetHighEnergyLimit(LivermoreHighEnergyLimit);
244 theComptonScattering->SetEmModel(theLivermoreComptonModel, 1);
245 ph->RegisterProcess(theComptonScattering, particle);
246
247 // gamma conversion - define low-energy model
248 G4GammaConversion* theGammaConversion = new G4GammaConversion();
249 G4VEmModel* theLivermoreGammaConversionModel =
251 theGammaConversion->SetEmModel(theLivermoreGammaConversionModel, 1);
252 ph->RegisterProcess(theGammaConversion, particle);
253
254 // default Rayleigh scattering is Livermore
255 G4RayleighScattering* theRayleigh = new G4RayleighScattering();
256 ph->RegisterProcess(theRayleigh, particle);
257
258 } else if (particleName == "e-") {
259
260 // multiple scattering
265 msc1->SetHighEnergyLimit(highEnergyLimit);
266 msc2->SetLowEnergyLimit(highEnergyLimit);
267 msc->AddEmModel(0, msc1);
268 msc->AddEmModel(0, msc2);
269
272 ss->SetEmModel(ssm, 1);
273 ss->SetMinKinEnergy(highEnergyLimit);
274 ssm->SetLowEnergyLimit(highEnergyLimit);
275 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
276 ph->RegisterProcess(msc, particle);
277 ph->RegisterProcess(ss, particle);
278
279 // Ionisation - Livermore should be used only for low energies
280 G4eIonisation* eIoni = new G4eIonisation();
281 G4LivermoreIonisationModel* theIoniLivermore = new
283 theIoniLivermore->SetHighEnergyLimit(0.1*MeV);
284 eIoni->AddEmModel(0, theIoniLivermore, new G4UniversalFluctuation() );
285 eIoni->SetStepFunction(0.2, 100*um); //
286 ph->RegisterProcess(eIoni, particle);
287
288 // Bremsstrahlung
290 G4VEmModel* theBremLivermore = new G4LivermoreBremsstrahlungModel();
291 theBremLivermore->SetHighEnergyLimit(1*GeV);
292 //theBremLivermore->SetAngularDistribution(new G4Generator2BS());
293 eBrem->SetEmModel(theBremLivermore,1);
294
295 ph->RegisterProcess(eBrem, particle);
296
297 } else if (particleName == "e+") {
298
299 // Identical to G4EmStandardPhysics_option3
300
301 // multiple scattering
306 msc1->SetHighEnergyLimit(highEnergyLimit);
307 msc2->SetLowEnergyLimit(highEnergyLimit);
308 msc->AddEmModel(0, msc1);
309 msc->AddEmModel(0, msc2);
310
313 ss->SetEmModel(ssm, 1);
314 ss->SetMinKinEnergy(highEnergyLimit);
315 ssm->SetLowEnergyLimit(highEnergyLimit);
316 ssm->SetActivationLowEnergyLimit(highEnergyLimit);
317
318 G4eIonisation* eIoni = new G4eIonisation();
319 eIoni->SetStepFunction(0.2, 100*um);
320
321 ph->RegisterProcess(msc, particle);
322 ph->RegisterProcess(eIoni, particle);
323 ph->RegisterProcess(new G4eBremsstrahlung(), particle);
324 ph->RegisterProcess(new G4eplusAnnihilation(), particle);
325 ph->RegisterProcess(ss, particle);
326
327 } else if (particleName == "mu+" ||
328 particleName == "mu-" ) {
329
330 G4MuIonisation* muIoni = new G4MuIonisation();
331 muIoni->SetStepFunction(0.2, 50*um);
332
333 ph->RegisterProcess(mumsc, particle);
334 ph->RegisterProcess(muIoni, particle);
335 ph->RegisterProcess(mub, particle);
336 ph->RegisterProcess(mup, particle);
337 ph->RegisterProcess(new G4CoulombScattering(), particle);
338
339 } else if (particleName == "alpha" ||
340 particleName == "He3" ) {
341
342 // Identical to G4EmStandardPhysics_option3
343
345 G4ionIonisation* ionIoni = new G4ionIonisation();
346 ionIoni->SetStepFunction(0.1, 10*um);
347
348 ph->RegisterProcess(msc, particle);
349 ph->RegisterProcess(ionIoni, particle);
350 ph->RegisterProcess(ionnuc, particle);
351
352 } else if (particleName == "GenericIon") {
353
354 // Identical to G4EmStandardPhysics_option3
355
356 G4ionIonisation* ionIoni = new G4ionIonisation();
358 ionIoni->SetStepFunction(0.1, 1*um);
359
360 ph->RegisterProcess(hmsc, particle);
361 ph->RegisterProcess(ionIoni, particle);
362 ph->RegisterProcess(ionnuc, particle);
363
364 } else if (particleName == "pi+" ||
365 particleName == "pi-" ) {
366
367 //G4hMultipleScattering* pimsc = new G4hMultipleScattering();
368 G4hIonisation* hIoni = new G4hIonisation();
369 hIoni->SetStepFunction(0.2, 50*um);
370
371 ph->RegisterProcess(pimsc, particle);
372 ph->RegisterProcess(hIoni, particle);
373 ph->RegisterProcess(pib, particle);
374 ph->RegisterProcess(pip, particle);
375
376 } else if (particleName == "kaon+" ||
377 particleName == "kaon-" ) {
378
379 //G4hMultipleScattering* kmsc = new G4hMultipleScattering();
380 G4hIonisation* hIoni = new G4hIonisation();
381 hIoni->SetStepFunction(0.2, 50*um);
382
383 ph->RegisterProcess(kmsc, particle);
384 ph->RegisterProcess(hIoni, particle);
385 ph->RegisterProcess(kb, particle);
386 ph->RegisterProcess(kp, particle);
387
388 } else if (particleName == "proton" ||
389 particleName == "anti_proton") {
390
391 //G4hMultipleScattering* pmsc = new G4hMultipleScattering();
392 G4hIonisation* hIoni = new G4hIonisation();
393 hIoni->SetStepFunction(0.2, 50*um);
394
395 ph->RegisterProcess(pmsc, particle);
396 ph->RegisterProcess(hIoni, particle);
397 ph->RegisterProcess(pb, particle);
398 ph->RegisterProcess(pp, particle);
399 ph->RegisterProcess(pnuc, particle);
400
401 } else if (particleName == "B+" ||
402 particleName == "B-" ||
403 particleName == "D+" ||
404 particleName == "D-" ||
405 particleName == "Ds+" ||
406 particleName == "Ds-" ||
407 particleName == "anti_He3" ||
408 particleName == "anti_alpha" ||
409 particleName == "anti_deuteron" ||
410 particleName == "anti_lambda_c+" ||
411 particleName == "anti_omega-" ||
412 particleName == "anti_sigma_c+" ||
413 particleName == "anti_sigma_c++" ||
414 particleName == "anti_sigma+" ||
415 particleName == "anti_sigma-" ||
416 particleName == "anti_triton" ||
417 particleName == "anti_xi_c+" ||
418 particleName == "anti_xi-" ||
419 particleName == "deuteron" ||
420 particleName == "lambda_c+" ||
421 particleName == "omega-" ||
422 particleName == "sigma_c+" ||
423 particleName == "sigma_c++" ||
424 particleName == "sigma+" ||
425 particleName == "sigma-" ||
426 particleName == "tau+" ||
427 particleName == "tau-" ||
428 particleName == "triton" ||
429 particleName == "xi_c+" ||
430 particleName == "xi-" ) {
431
432 // Identical to G4EmStandardPhysics_option3
433
434 ph->RegisterProcess(hmsc, particle);
435 ph->RegisterProcess(new G4hIonisation(), particle);
436 ph->RegisterProcess(pnuc, particle);
437 }
438 }
439
440 // Em options
441 //
443 opt.SetVerbose(verbose);
444
445 // Multiple Coulomb scattering
446 //
447 opt.SetPolarAngleLimit(CLHEP::pi);
448
449 // Physics tables
450 //
451
452 opt.SetMinEnergy(10*eV);
453 opt.SetMaxEnergy(10*TeV);
454 opt.SetDEDXBinning(240);
455 opt.SetLambdaBinning(240);
456
457 // Nuclear stopping
458 pnuc->SetMaxKinEnergy(MeV);
459
460 // Ionization
461 //
462 //opt.SetSubCutoff(true);
463
464 // Deexcitation
465 //
468 de->SetFluo(true);
469}
@ 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 SetActivationLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:606
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:592
void SetMinKinEnergy(G4double e)
void SetEmModel(G4VEmModel *, G4int index=1)
void SetMaxKinEnergy(G4double e)
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: