Penelope low-energy electromagnetic models
From version 9.2 of Geant4, the Penelope low-energy electromagnetic models for e± and gamma-rays have been migrated to the new design (the same as for "Standard EM processes"). From version Geant4 9.4 a complete new set of Penelope models has been made available, which is based on the version 2008 of Penelope, rather than version 2001.
All Penelope models are compliant with the multi-thread approach employed in Geant4 from version 10.0.
- The model approach and the available Penelope models
- How to use the new models in Geant4 user applications
- Additional controls/settings for the models
- How to access/retrieve physics information from the models
- How to simulate atomic deexcitation?
- Production cuts
- Validation and Papers
The model approach and the available Penelope model
Physics process |
Particle(s) |
Old Penelope Process | New Penelope Model |
Rayleigh scattering | gamma | G4PenelopeRayleigh | G4PenelopeRayleighModel |
Compton scattering | gamma | G4PenelopeCompton | G4PenelopeComptonModel |
Photo-electric effect | gamma | G4PenelopePhotoElectric | G4PenelopePhotoElectricModel |
Pair production | gamma | G4PenelopeGammaConversion | G4PenelopeGammaConversionModel |
Ionisation | e± | G4PenelopeIonisation | G4PenelopeIonisationModel |
Bremsstrahlung | e± | G4PenelopeBremsstrahlung | G4PenelopeBremsstrahlungModel |
Positron annihilation | e+ | G4PenelopeAnnihilation | G4PenelopeAnnihilationModel |
In the old approach there were multiple versions of the same process, one in the standard EM package and one/two in the LowEnergy EM package (e.g. for Compton Scattering: G4LowEnergyCompton, G4PenelopeCompton and G4ComptonScattering).
In the new approach there is only one process (e.g. G4ComptonScattering) and multiple models that can be registered to the process (same approch as hadronic physics). For instance, the models available for the Compton scattering (to be registered to G4ComptonScattering) are: G4KleinNishinaCompton (default), G4PenelopeComponModel(Penelope) and G4LivermoreComptonModel (former LowEnergy model).
How to use the new models in Geant4 user applications
The most direct way to implement the Penelope physics in a user application is to register the ready-for-the-use constructor G4EmPenelopePhysics provided by Geant4 into the application physics list. The constructor can be directly registered to the modular physics list of the user application as
#include "G4EmPenelopePhysics.hh";
G4VPhysicsConstructor* emPhysicsList = new G4EmPenelopePhysics();
emPhysicsList->ConstructProcess();
or
#include "G4EmPenelopePhysics.hh";
RegisterPhysics( new G4EmPenelopePhysics() );
To get the full flexibility about the energy window and the types of models to be registered for each physics process, the application physics list must be written ad-hoc. The following code can be used in the user physics list, to register the individual Penelope models to the EM processes: first instantiate the process and then register the model to the process.
In this example, the Penelope model is used up to 1 GeV; above 1 GeV the default model is considered (which is G4KleinNishinaCompton for Compton Scattering and G4MollerBhabhaModel for electron ionisation).
#include "G4PenelopeComptonModel.hh";
#include "G4ComptonScattering.hh";
G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
theComptonScattering->SetEmModel(new G4PenelopeComptonModel());
theComptonScattering->SetHighEnergy(1.*GeV);
pmanager->AddDiscreteProcess(theComptonScattering);
#include "G4PenelopeIonisationModel.hh";
#include "G4eIonisation.hh";
G4eIonisation* theIonisation = new G4eIonisation();
theIonisation->SetEmModel(new G4PenelopeIonisationModel());
heIonisation->SetHighEnergy(1.*GeV);
pmanager->AddProcess(theIonisation,-1,1,1);
The processes to which the Penelope models have to be registered are summarized in the table below:
Physics process |
Particle(s) |
Process to which register the model | Penelope Model |
Rayleigh scattering | gamma | G4RayleighScattering | G4PenelopeRayleighModel |
Compton scattering | gamma | G4ComptonScattering | G4PenelopeComptonModel |
Photo-electric effect | gamma | G4PhotoElectricEffect | G4PenelopePhotoElectricModel |
Pair production | gamma | G4GammaConversion | G4PenelopeGammaConversionModel |
Ionisation | e± | G4eIonisation | G4PenelopeIonisationModel |
Bremsstrahlung | e± | G4eBremsstrahlung | G4PenelopeBremsstrahlungModel |
Positron annihilation | e+ | G4eplusAnnihilation | G4PenelopeAnnihilationModel |
The ordering of the continuous processes for electrons and positrons is the same as for Standard and Livermore models, namely:
electrons:
pmanager->AddProcess(new G4eMultipleScattering,-1, 1, 1);
pmanager->AddProcess(new G4eIonisation, -1, 2, 2);
pmanager->AddProcess(new G4eBremsstrahlung(), -1,-3, 3);
positrons:
pmanager->AddProcess(new G4eMultipleScattering,-1, 1, 1);
pmanager->AddProcess(new G4eIonisation, -1, 2, 2);
pmanager->AddProcess(new G4eBremsstrahlung, -1,-3, 3);
pmanager->AddProcess(new G4eplusAnnihilation, 0,-1, 4);
Notice that at the moment, the G4Penelope package does not include an alternative model for G4eMultipleScattering.
Additional controls/settings for the models
All Penelope models have a method called SetVerbosityLevel() to set the verbosity level (default verbosity is 0).
The verbosity scales goes from 0 to 4:
0 = no verbosity
1 = warning for energy non-conservation
2 = details of energy budget
3 = calculation of cross sections/stopping powers, file openings, sampling of atoms
4 = entering in methods
How to access/retrieve physics information from the models
The new model approach makes the physics content more transparent to users. It also makes possible to use all general tools (e.g. G4EmCalculator) that are available for the "Standard" electromagnetic models.
For instance, this is the code to retrieve directly from the model the mean free path and the stopping power from G4PenelopeIonisationModel for electrons in a given material.
The code below gives: (1) the ionisation mean free path in the material aMaterial for production of delta rays above energy=energyCut for electrons of energy initEnergy; (2) the stopping power due to soft ionisation (below energyCut).
#include "G4PenelopeIonisationModel.hh"
#include "G4Electron.hh"
G4PenelopeIonisationModel* ioni = new G4PenelopeIonisationModel();
G4DataVector dummy;
ioni->Initialise(G4Electron::Definition(),dummy); //initialise the model
ioni->SetVerbosityLevel(1); //optional verbosity level
//
G4double inverseMFP = ioni->CrossSectionPerVolume(aMaterial,
G4Electron::Definition(), initEnergy,energyCut);
G4cout << "Mean free path in " << aMaterial->GetName() <<
" at " << initEnergy/keV <<
" keV for energy > " << energyCut/keV << " keV = " <<
(1./inverseMFP)/mm << " mm" << G4endl;
//
G4double sPower = ioni->ComputeDEDXPerVolume(aMaterial,
G4Electron::Definition(), initEnergy,energyCut);
G4cout << "Stopping power path in " << aMaterial->GetName() <<
" at " << initEnergy/keV <<
" keV for energy < " << energyCut/keV << " keV = " <<
sPower/(keV/mm) << "keV/mm" << G4endl;
//
For more examples, see also the extended example electromagnetic/TestEm0
How to simulate atomic deexcitation?
Follow this page.
Production cuts
Remember that production cuts for secondaries can be specified as range cuts, which at initialisation time are converted into energy threshold for secondary gamma, electron, positron and proton production.
1) Range cut value is set to 0.7 mm in Geant4 reference Physics Lists. This value can be specified in the optional SetCuts() method of your Physics list or via UI command, for eg. to set a range cut of 10 micrometers, one can use:
/run/setCut 0.01 mm
or for a given particle type (for e.g. electron)
/run/setCutForAGivenParticle e- 0.01 mm
2) Since not all Geant4 models are able to work with very low production thresholds, an energy threshold limit is used, its default value is set to 990 eV. You can change this value (for eg. to 250 eV) by using the UI command:
/cuts/setLowEdge 250 eV
or alternatively directly in your Physics list in the optional SetCuts() method with:
G4ProductionCutsTable::GetProductionCutsTable()->SetEnergyRange(250*eV, 1*GeV);
Validation and papers
Bremsstrahlung: please refer to Nucl.Instrum.Meth. B350 (2015) 41-48 for a quantitative comparison of bremsstrahlung models at low energy (below a few MeV). All bremsstrahlung models are more or less equivalent in the simulation of absolute double-differential energy/angular distributions in thick targets. Still, the Penelope model is doing globally slightly better than Livermore and Opt3.