The Geant4-DNA project
Contents of this page:
Purpose
Geant4 is being extended with very low energy processes for the
modelling of biological damages induced by ionising radiation. Such developments are on-going in the framework of the
Geant4-DNA project, initiated by in 2000
P. Nieminen (
European Space Agency/ESTEC
) in order to model physics, chemistry and biology processes for the simulation of early biological effects of ionizing radiation.
The
first phase of the project was undertaken by INFN Genova in Italy, delivering Work Package reports and a User Requirement Document. The
second phase is currently on-going since 2004. It is now coordinated by
S. Incerti and it is a
full activity of the Low Energy Electromagnetic Physics working group of the Geant4 collaboration.
The
rules and regulations of the Geant4-DNA collaboration as well as the list of collaborators are available from this
page.
Expected developments include :
- Physics processes in liquid water and other biological materials
- Chemistry and physico-chemistry processes
- Molecular geometries
- Quantification of damages (single-strand, double-strand breaks, ...)
Publications
The Geant4-DNA Physics processes are described in details in the following publications :
- Electron and proton elastic scattering in water vapour, C. Champion, S. Incerti, H. Tran, Z. El Bitar, Nucl. Instrum. and Meth. B (2011) in press (link
)
- Combination of electromagnetic Physics processes for microdosimetry in liquid water with the Geant4 Monte Carlo simulation toolkit, V.N. Ivanchenko, S. Incerti, Z. Francis, H.N. Tran, M. Karamitros, M.A. Bernal, C. Champion and P. Guèye, Nucl. Instrum. and Meth. B (2011) in press (link
)
- Modeling radiation chemistry in the Geant4 toolkit, M. Karamitros, A. Mantero, S. Incerti et al., in press
- Recent improvements in Geant4 electromagnetic physics models and interfaces, V. Ivantchenko et al., in press
- The invariance of the total direct DNA strand break yield, M. A. Bernal, C. E. de Almeida, C. S. Sampaio, S. Incerti, C. Champion, P. Nieminen, Med. Phys., in press (link
)
- Stopping power and ranges of electrons, protons and alpha particles in liquid water using the Geant4-DNA package, Z. Francis, S. Incerti, M. Karamitros, H.N. Tran, C. Villagrasa, Nucl. Instrum. and Meth. B (2011) in press (link
)
- Molecular scale track structure simulations in liquid water using the Geant4-DNA Monte Carlo processes, Z. Francis, S. Incerti, R. Capra, B. Mascialino, G. Montarou, V. Stepan, C. Villagrasa, Appl. Radiat. Isot. 69 (2011) 220-226 (link
)
- Physical models implemented in the Geant4-DNA extension of the Geant4 toolkit for calculating initial radiation damage at the molecular level, C. Villagrasa, Z. Francis, S. Incerti, Rad. Prot. Dos. 133, 1 (2010) 1-5 (link
)
- Comparison of Geant4 very low energy cross section models with experimental data in water, S. Incerti, A. Ivanchenko, M. Karamitros, A. Mantero, P. Moretto, H. N. Tran, B. Mascialino, C. Champion, V. N. Ivanchenko, M. A. Bernal, Z. Francis, C. Villagrasa, G. Baldacchino, P. Guèye, R. Capra, P. Nieminen, C. Zacharatou, Med. Phys. 37 (2010) 4692-4708 (link
)
- A free-parameter theoretical model for describing the electron elastic scattering in water in the Geant4 toolkit, C. Champion, S. Incerti, H. Aouchiche, D. Oubaziz, Rad. Phys. Chem. 78 (2009) 745-750 (link
)
- The Geant4-DNA project, S. Incerti, G. Baldacchino, M. Bernal, R. Capra, C. Champion, Z. Francis, S. Guatelli, P. Guèye, A. Mantero, B. Mascialino, P. Moretto, P. Nieminen, A. Rosenfeld, C. Villagrasa and C. Zacharatou, Int. J. Model. Simul. Sci. Comput. 1 (2010) 157–178 (link
)
Our
full list of publications including activities at the Physics-Biology frontier are available from
this page
.
Design of Physics process and model classes
The Geant4-DNA design has been adopted by the Standard Electromagnetic Physics working group in 2006. A single physical process can handle several models which are registered to the process.
- All Geant4-DNA process classes derive from the G4VEmProcess abstract base class, which derives from the G4VDiscreteProcess class, since all Geant4-DNA processes are discrete processes. They use the naming format: G4DNAXXX where XXX is the process name (for eg. G4DNAExcitation for excitation). This is the name that will be displayed when using the command /tracking/verbose 1. Several methods are implemented :
- G4DNAXXX::InitialiseProcess(...) : process initialisation
- G4DNAXXX::PrintInfo(...) : verbose
- it is also possible to customize the process name : when a process is added to the process manager in the Physics list (see example below), one can use a user defined name, for eg. "myName" :
pmanager->AddDiscreteProcess(new G4DNAExcitation("myName"));
- All Geant4-DNA model classes derive from the G4VEmModel abstract base class. They use the naming format: G4DNAYYYModel, where YYY is the model name (for eg. G4DNAEmfietzoglouExcitationModel for excitation). Several methods are implemented :
- G4DNAYYYModel::Initialise(...) : model initialisation
- G4DNAYYYModel::CrossSectionPerVolume(...) : computation of model total cross section multiplied by material density (equivalent to inverse mean free path)
- G4DNAYYYModel::SampleSecondaries(...) : computation of model final state
Available Physics processes and models
We indicate in the table below the
Geant4-DNA physics processes, models and particles available in the Geant4 toolkit from June 2009, as well as their
energy range of applicability. Each process can handle several alternative or complementary (in energy) models.
Material |
Physics process |
Process class |
Model class |
Low enery applicability limit of model class |
High enery applicability limit of model class |
Model class kills incident particle below this energy |
Data tables for interpolated models (in $G4LEDATA/dna) |
|
e- |
Liquid water |
Elastic scattering |
G4DNAElastic |
G4DNAChampionElasticModel (default model in the G4EmDNAPhysics builder - see Physics builders) |
0 eV |
1 MeV |
4 eV - can be decreased by user (see below) but not validated |
sigma_elastic_e_champion.dat, sigmadiff_elastic_e_champion.dat |
Liquid water |
Elastic scattering |
G4DNAElastic |
G4DNAScreenedRutherfordElasticModel (alternative model) |
0 eV |
1 MeV |
9 eV - can be decreased by user (see below) but not validated |
n/a (analytical model) |
Liquid water |
Excitation |
G4DNAExcitation |
G4DNABornExcitationModel |
9 eV |
1MeV |
n/a |
sigma_excitation_e_born.dat |
Liquid water |
Ionisation |
G4DNAIonisation |
G4DNABornIonisationModel |
11 eV |
1 MeV |
n/a |
sigma_ionisation_e_born.dat, sigmadiff_ionisation_e_born.dat |
Liquid water |
Vib. excitation |
G4DNAVibExcitation |
G4DNASancheExcitationModel |
2 eV - can be decreased by user (see below) but not validated |
100 eV |
n/a |
sigma_excitationvib_e_sanche.dat |
Liquid water |
Attachment |
G4DNAAttachment |
G4DNAMeltonAttachmentModel |
4 eV |
13 eV |
n/a |
sigma_attachment_e_melton.dat |
proton |
Liquid water |
Excitation |
G4DNAExcitation |
G4DNAMillerGreenExcitationModel |
10 eV |
500 keV |
n/a |
n/a (analytical model) |
Liquid water |
Excitation |
G4DNAExcitation |
G4DNABornExcitationModel |
500 keV |
100 MeV |
n/a |
sigma_excitation_p_born.dat |
Liquid water |
Ionisation |
G4DNAIonisation |
G4DNARuddIonisationModel |
0 eV |
500 keV |
100 eV |
sigma_ionisation_p_rudd.dat |
Liquid water |
Ionisation |
G4DNAIonisation |
G4DNABornIonisationModel |
500 keV |
100 MeV |
n/a |
sigma_ionisation_p_born.dat, sigmadiff_ionisation_p_born.dat |
Liquid water |
Charge decrease |
G4DNAChargeDecrease |
G4DNADingfelderChargeDecreaseModel |
100 eV |
10 MeV |
n/a |
n/a (analytical model) |
hydrogen (neutral atom) |
Liquid water |
Excitation |
G4DNAExcitation |
G4DNAMillerGreenExcitationModel |
10 eV |
500 keV |
n/a |
n/a (analytical model) |
Liquid water |
Ionisation |
G4DNAIonisation |
G4DNARuddIonisationModel |
0 eV |
100 MeV |
100 eV |
sigma_ionisation_h_rudd.dat |
Liquid water |
Charge increase |
G4DNAChargeIncrease |
G4DNADingfelderChargeIncreaseModel |
100 eV |
10 MeV |
n/a |
n/a (analytical model) |
alpha (neutral Helium atom ionised twice) |
Liquid water |
Excitation |
G4DNAExcitation |
G4DNAMillerGreenExcitationModel |
1 keV |
400 MeV |
n/a |
n/a (analytical model) |
Liquid water |
Ionisation |
G4DNAIonisation |
G4DNARuddIonisationModel |
0 keV |
400 MeV |
1 keV |
sigma_ionisation_alphaplusplus_rudd.dat |
Liquid water |
Charge decrease |
G4DNAChargeDecrease |
G4DNADingfelderChargeDecreaseModel |
1 keV |
10 MeV |
n/a |
n/a (analytical model) |
alpha+ (neutral Helium atom ionised once) |
Liquid water |
Excitation |
G4DNAExcitation |
G4DNAMillerGreenExcitationModel |
1 keV |
400 MeV |
n/a |
n/a (analytical model) |
Liquid water |
Ionisation |
G4DNAIonisation |
G4DNARuddIonisationModel |
0 keV |
400 MeV |
1 keV |
sigma_ionisation_alphaplus_rudd.dat |
Liquid water |
Charge decrease |
G4DNAChargeDecrease |
G4DNADingfelderChargeDecreaseModel |
1 keV |
10 MeV |
n/a |
n/a (analytical model) |
Liquid water |
Charge increase |
G4DNAChargeIncrease |
G4DNADingfelderChargeIncreaseModel |
1 keV |
10 MeV |
n/a |
n/a (analytical model) |
helium (neutral Helium atom) |
Liquid water |
Excitation |
G4DNAExcitation |
G4DNAMillerGreenExcitationModel |
1 keV |
400 MeV |
n/a |
n/a (analytical model) |
Liquid water |
Ionisation |
G4DNAIonisation |
G4DNARuddIonisationModel |
0 keV |
400 MeV |
1 keV |
sigma_ionisation_he_rudd.dat |
Liquid water |
Charge increase |
G4DNAChargeIncrease |
G4DNADingfelderChargeIncreaseModel |
1 keV |
10 MeV |
n/a |
n/a (analytical model) |
C, N, O, Fe |
Liquid water |
Ionisation |
G4DNAIonisation |
G4DNARuddIonisationExtendedModel (preliminary) - must be registered by the user to G4DNAIonisation |
0.5 MeV/u |
1e6 MeV/u |
0.5 MeV/u |
sigma_ionisation_c_rudd.dat, sigma_ionisation_n_rudd.dat, sigma_ionisation_o_rudd.dat, sigma_ionisation_fe_rudd.dat |
Example Physics list
We indicate below a
Physics list example for the use of
Geant4-DNA physics processes and models only.
...
#include "G4DNAGenericIonsManager.hh"
...
void PhysicsList::ConstructBaryons()
{
...
// Geant4-DNA particles
G4GenericIon::GenericIonDefinition();
G4DNAGenericIonsManager * genericIonsManager;
genericIonsManager=G4DNAGenericIonsManager::Instance();
genericIonsManager->GetIon("alpha++");
genericIonsManager->GetIon("alpha+");
genericIonsManager->GetIon("helium");
genericIonsManager->GetIon("hydrogen");
genericIonsManager->GetIon("carbon");
genericIonsManager->GetIon("nitrogen");
genericIonsManager->GetIon("oxygen");
genericIonsManager->GetIon("iron");
...
}
...
#include "G4DNAElastic.hh"
#include "G4DNAChampionElasticModel.hh"
#include "G4DNAScreenedRutherfordElasticModel.hh"
#include "G4DNAExcitation.hh"
#include "G4DNAMillerGreenExcitationModel.hh"
#include "G4DNABornExcitationModel.hh"
#include "G4DNAIonisation.hh"
#include "G4DNABornIonisationModel.hh"
#include "G4DNARuddIonisationModel.hh"
#include "G4DNAChargeDecrease.hh"
#include "G4DNADingfelderChargeDecreaseModel.hh"
#include "G4DNAChargeIncrease.hh"
#include "G4DNADingfelderChargeIncreaseModel.hh"
#include "G4DNAVibExcitation.hh"
#include "G4DNASancheExcitationModel.hh"
#include "G4DNAAttachment.hh"
#include "G4DNAMeltonAttachmentModel.hh"
...
void PhysicsList::ConstructEM()
{
theParticleIterator->reset();
while( (*theParticleIterator)() )
{
G4ParticleDefinition* particle = theParticleIterator->value();
G4ProcessManager* pmanager = particle->GetProcessManager();
G4String particleName = particle->GetParticleName();
if (particleName == "e-") {
G4DNAElastic* theDNAElasticProcess = new G4DNAElastic("e-_G4DNAElastic");
theDNAElasticProcess->SetModel(new G4DNAChampionElasticModel());
// or alternative model (faster but much less accurate at low energies)
// theDNAElasticProcess->SetModel(new G4DNAScreenedRutherfordElasticModel());
pmanager->AddDiscreteProcess(theDNAElasticProcess);
pmanager->AddDiscreteProcess(new G4DNAExcitation("e-_G4DNAExcitation"));
pmanager->AddDiscreteProcess(new G4DNAIonisation("e-_G4DNAIonisation"));
pmanager->AddDiscreteProcess(new G4DNAVibExcitation("e-_G4DNAVibExcitation"));
pmanager->AddDiscreteProcess(new G4DNAAttachment("e-_G4DNAAttachment"));
} else if ( particleName == "proton" ) {
pmanager->AddDiscreteProcess(new G4DNAExcitation("proton_G4DNAExcitation"));
pmanager->AddDiscreteProcess(new G4DNAIonisation("proton_G4DNAIonisation"));
pmanager->AddDiscreteProcess(new G4DNAChargeDecrease("proton_G4DNAChargeDecrease"));
} else if ( particleName == "hydrogen" ) {
pmanager->AddDiscreteProcess(new G4DNAExcitation("hydrogen_G4DNAExcitation"));
pmanager->AddDiscreteProcess(new G4DNAIonisation("hydrogen_G4DNAIonisation"));
pmanager->AddDiscreteProcess(new G4DNAChargeIncrease("hydrogen_G4DNAChargeIncrease"));
} else if ( particleName == "alpha" ) {
pmanager->AddDiscreteProcess(new G4DNAExcitation("alpha_G4DNAExcitation"));
pmanager->AddDiscreteProcess(new G4DNAIonisation("alpha_G4DNAIonisation"));
pmanager->AddDiscreteProcess(new G4DNAChargeDecrease("alpha_G4DNAChargeDecrease"));
} else if ( particleName == "alpha+" ) {
pmanager->AddDiscreteProcess(new G4DNAExcitation("alpha+_G4DNAExcitation"));
pmanager->AddDiscreteProcess(new G4DNAIonisation("alpha+_G4DNAIonisation"));
pmanager->AddDiscreteProcess(new G4DNAChargeDecrease("alpha+_G4DNAChargeDecrease"));
pmanager->AddDiscreteProcess(new G4DNAChargeIncrease("alpha+_G4DNAChargeIncrease"));
} else if ( particleName == "helium" ) {
pmanager->AddDiscreteProcess(new G4DNAExcitation("helium_G4DNAExcitation"));
pmanager->AddDiscreteProcess(new G4DNAIonisation("helium_G4DNAIonisation"));
pmanager->AddDiscreteProcess(new G4DNAChargeIncrease("helium_G4DNAChargeIncrease"));
}
} // Loop on particles
}
Physics builders
Instead of creating your own Physics list, we recommend the usage of Physics builders. The
G4EmDNAPhysics builder is available in
$G4INSTALL/source/physics_lists/builders.
This builder contains the same Physics list as above.
For example, the
dnaphysics advanced example located in
$G4INSTALL/examples/advanced shows how to use the Geant4-DNA Physics builder.
You may also look at the following
Physics list that explains how to use a G4EmDNAPhysics object.
How to access Physics ?
It is possible to retrieve Physics quantities using a
G4EmCalculator object. For example, in order to retrieve the
total cross section of a process with name procName, do as follows :
#include "G4EmCalculator.hh"
...
G4EmCalculator emCalculator;
G4double density = material->GetDensity();
G4double massSigma = emCalculator.ComputeCrossSectionPerVolume(energy,particle,procName,material)/density;
G4cout << G4BestUnit(massSigma, "Surface/Mass") << G4endl;
A good example is
TestEm14 located in
$G4INSTALL/examples/extended/electromagnetic, look in particular at the
RunAction.cc class
Note that in order to use the
ComputeCrossSectionPerVolume method from the
G4EmCalculator class for
alpha particles (G4Alpha), it is necessary to add the following condition to
G4EmCalculator::UpdateParticle in
G4EmCalculator.cc (located in $G4INSTALL/source/processes/electromagnetic/utils/src)
&& currentParticleName != "alpha"
How to kill particles below a given energy threshold for faster performance ?
In case you want to kill particles with energies below an energy threshold value, you may instantiate a
G4UserLimits object in your
DetectorConstruction class and define the process
G4UserSpecialCuts in your
PhysicsList class for the affected particles. All details are given in the
Geant4 User's Guide For Application Developers
.
For example, if you want to kill all electrons below 9 eV, you may use the following lines. All electron tracks below 9 eV will be killed and electrons will deposit locally their total energy :
- in your
DetectorConstruction class, if you choose to apply this limit to the World volume :
#include "G4UserLimits.hh"
...
logicWorld->SetUserLimits(new G4UserLimits(DBL_MAX,DBL_MAX,DBL_MAX,9*eV));
- in your
PhysicsList class :
#include "G4UserSpecialCuts.hh"
...
if (particleName == "e-")
{
...
pmanager->AddDiscreteProcess(new G4UserSpecialCuts());
...
}
How to decrease the low energy limit under which electron elastic scattering and vibrational excitation models kill electrons ?
Please do this in your own PhysicsList
at your own risk since these models are not experimentally validated at very low energies. You will receive
WARNINGS at initialisation.
// Elastic scattering (can be extended down to 0.025 eV)
G4DNAElastic* theDNAElasticProcess = new G4DNAElastic("e-_G4DNAElastic");
theDNAElasticProcess->SetModel(new G4DNAChampionElasticModel());
// or alternative model
// theDNAElasticProcess->SetModel(new G4DNAScreenedRutherfordElasticModel());
((G4DNAChampionElasticModel*)(theDNAElasticProcess->Model()))->SetKillBelowThreshold(0.025*eV);
// or alternatively
// ((G4DNAScreenedRutherfordElasticModel*)(theDNAElasticProcess->Model()))->SetKillBelowThreshold(0.025*eV);
pmanager->AddDiscreteProcess(theDNAElasticProcess);
// Vibrational excitation (can be extended down to 0.025 eV)
pmanager->AddDiscreteProcess(new G4DNAVibExcitation("e-_G4DNAVibExcitation"));
G4DNAVibExcitation* theDNAVibProcess = new G4DNAVibExcitation("e-_G4DNAVib");
theDNAVibProcess->SetModel(new G4DNASancheExcitationModel() );
((G4DNASancheExcitationModel*)(theDNAVibProcess->Model()))->ExtendLowEnergyLimit(0.025*eV);
pmanager->AddDiscreteProcess(theDNAVibProcess);
How to activate Auger production ?
Put the following in your macro file,
provided that you use the G4DNAEmPhysics Physics builder.
/run/initialize
/process/em/auger true
Alternatively, consult this
link.
Material definition and density change
Geant4-DNA Physics models can only be used in liquid water defined as the NIST material
G4_WATER. To do so, you may use the following lines in your detector construction class :
// Water is defined from NIST material database
G4NistManager * man = G4NistManager::Instance();
G4Material * H2O = man->FindOrBuildMaterial("G4_WATER");
// Default materials in setup.
waterMaterial = H2O;
In case you
need to change the density of water, you may look into the
dnaphysics advanced example DetectorConstruction class:
// If one wishes to test other density value for water material, one should use instead:
G4Material * H2O = man->BuildMaterialWithNewDensity("G4_WATER_MODIFIED","G4_WATER",1000*g/cm/cm/cm);
// Note: any string for "G4_WATER_MODIFIED" parameter is accepted and "G4_WATER" parameter should not be changed
Geometrical tolerances
When using Geant4-DNA processes, the geometrical tolerances of your setup must be adapted to the size of your setup when working with nanometer size dimensions. For this, in your DetectorConstruction class construct() method, you may use for eg. the following lines:
#include "G4GeometryManager.hh"
#include "G4GeometryTolerance.hh"
...
G4double fWorldLength = 100. *nm;
G4GeometryManager::GetInstance()->SetWorldMaximumExtent(fWorldLength);
G4cout << "Surface tolerance = " << G4GeometryTolerance::GetInstance()->GetSurfaceTolerance()/nm
<< " nm - Radial tolerance = " << G4GeometryTolerance::GetInstance()->GetRadialTolerance()/nm << " nm" <<G4endl;
Combination with other Physics processes for multi-scale simulations
The new
microdosimetry advanced example (
link
) shows how to use these Geant4-DNA processes in combination with Standard Electromagnetic Physics processes.
Examples ready-to-use
Some examples are already available in Geant4 for you to test Geant4-DNA processes and models.
Example code name |
Purpose |
Where to find it ? |
Comments |
dnaphysics |
usage of Geant4-DNA physics builder |
$G4INSTALL/examples/advanced |
from Geant4 9.5 BETA |
microdosimetry |
combination of Geant4-DNA processes with other EM processes |
$G4INSTALL/examples/advanced |
from Geant4 9.5 BETA |
TestEm7 |
usage of Physics builders (see PhysicsList class) |
$G4INSTALL/examples/extended |
from Geant4 9.4 |
TestEm12 |
dosimetry in spherical shells with Geant4-DNA |
$G4INSTALL/examples/extended |
from Geant4 9.5 BETA |
TestEm14 |
extraction of cross section (see RunAction class) |
$G4INSTALL/examples/extended |
from Geant4 9.4 |
The corresponding advanced examples are documented
here
.
If you do not wish to install yourself Geant4, you may freely use our Geant4 Virtual Machine for VmWare and VirtualBox, downloadable from this link
.
Main Low Energy WG page