|
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").
Notice: the new Penelope models are released as
"beta" software. Deeper test and validation is ongoing.
Feedback and comments from users are very welcome to improve the new code to production quality and to
report for possible bugs.
The original Penelope processes are momentarily kept for
backwards compatibility, but they will be removed in next Geant4 releases.
At the moment, the physics results by the original G4Penelope processes and the newly migrated
models are exactly the same.
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
The model approach and the available Penelope models
| 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 following code can be used in the physics list, to
register the Penelope models to the EM processes:
first instantiate the process and then register the model to the process.
#include "G4PenelopeComptonModel.hh";
#include "G4ComptonScattering.hh";
G4ComptonScattering* theComptonScattering = new G4ComptonScattering();
theComptonScattering->SetModel(new G4PenelopeComptonModel());
pmanager->AddDiscreteProcess(theComptonScattering);
#include "G4PenelopeIonisationModel.hh";
#include "G4eIonisation.hh";
G4eIonisation* theIonisation = new G4eIonisation();
theIonisation->SetModelEm(new G4PenelopeIonisationModel());
pmanager->AddProcess(theIonisation,-1,1,1);
Notice: for discrete processes (Compton and Rayleigh
scattering, pair production, photo-electric effect and positron
annihilation), the method to invoke is
SetModel()
For continuous-discrete processes (bremsstrahlung and ionisation), the
method is
SetEmModel().
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 |
g |
G4RayleighScattering |
G4PenelopeRayleighModel |
|
Compton scattering |
g |
G4ComptonScattering |
G4PenelopeComptonModel |
|
Photo-electric effect |
g |
G4PhotoElectricEffect |
G4PenelopePhotoElectricModel |
|
Pair production |
g |
G4GammaConversion |
G4PenelopeGammaConversionModel |
|
Ionisation |
e± |
G4eIonisation |
G4PenelopeIonisationModel |
|
Bremsstrahlung |
e± |
G4eBremsstrahlung |
G4PenelopeBremsstrahlungModel |
|
Positron annihilation |
e+ |
G4eplusAnnihilation |
G4PenelopeAnnihilationModel |
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
Penelope models which are interfaced to the G4AtomicDeexcitation (ionisation,
Compton and photo-electric) have a method called
SetUseAtomicDeexcitation()
to switch on/off the fluorescence
production. Default is true,
namely fluorescence is explicitely produced. If it set to
false, the energy for
fluorescence products is deposited locally.
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
|