Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
|
The G4VDNAModel class. More...
#include <G4VDNAModel.hh>
Public Member Functions | |
G4VDNAModel (const G4String &nam, const G4String &applyToMaterial) | |
G4VDNAModel Constructeur of the G4VDNAModel class. | |
virtual | ~G4VDNAModel () |
~G4VDNAModel | |
virtual void | Initialise (const G4ParticleDefinition *particle, const G4DataVector &cuts, G4ParticleChangeForGamma *fpChangeForGamme=nullptr)=0 |
Initialise Each model must implement an Initialize method. | |
virtual G4double | CrossSectionPerVolume (const G4Material *material, const G4String &materialName, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)=0 |
CrossSectionPerVolume Every model must implement its own CrossSectionPerVolume method. It is used by the process to determine the step path and must return a cross section times a number of molecules per volume unit. | |
virtual void | SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4String &materialName, const G4DynamicParticle *, G4ParticleChangeForGamma *particleChangeForGamma, G4double tmin=0, G4double tmax=DBL_MAX)=0 |
SampleSecondaries Each model must implement SampleSecondaries to decide if a particle will be created after the ModelInterface or if any charateristic of the incident particle will change. | |
G4bool | IsMaterialDefine (const G4String &materialName) |
IsMaterialDefine Check if the given material is defined in the simulation. | |
G4bool | IsMaterialExistingInModel (const G4String &materialName) |
IsMaterialExistingInModel Check if the given material is defined in the current model class. | |
G4bool | IsParticleExistingInModelForMaterial (const G4String &particleName, const G4String &materialName) |
IsParticleExistingInModelForMaterial To check two things: 1- is the material existing in model ? 2- if yes, is the particle defined for that material ? | |
G4String | GetName () |
GetName. | |
G4double | GetHighELimit (const G4String &material, const G4String &particle) |
GetHighEnergyLimit. | |
G4double | GetLowELimit (const G4String &material, const G4String &particle) |
GetLowEnergyLimit. | |
void | SetHighELimit (const G4String &material, const G4String &particle, G4double lim) |
SetHighEnergyLimit. | |
void | SetLowELimit (const G4String &material, const G4String &particle, G4double lim) |
SetLowEnergyLimit. | |
Protected Types | |
typedef std::map< G4String, std::map< G4String, G4DNACrossSectionDataSet *, std::less< G4String > > > | TableMapData |
typedef std::map< G4String, std::map< G4String, G4double > > | RatioMapData |
typedef std::map< G4String, G4double >::const_iterator | ItCompoMapData |
Protected Member Functions | |
TableMapData * | GetTableData () |
GetTableData. | |
std::vector< G4String > | BuildApplyToMatVect (const G4String &materials) |
BuildApplyToMatVect Build the material name vector which is used to know the materials the user want to include in the model. | |
void | ReadAndSaveCSFile (const G4String &materialName, const G4String &particleName, const G4String &file, G4double scaleFactor) |
ReadAndSaveCSFile Read and save a "simple" cross section file : use of G4DNACrossSectionDataSet->loadData() | |
G4int | RandomSelectShell (G4double k, const G4String &particle, const G4String &materialName) |
RandomSelectShell Method to randomely select a shell from the data table uploaded. The size of the table (number of columns) is used to determine the total number of possible shells. | |
void | AddCrossSectionData (G4String materialName, G4String particleName, G4String fileCS, G4String fileDiffCS, G4double scaleFactor) |
AddCrossSectionData Method used during the initialization of the model class to add a new material. It adds a material to the model and fills vectors with informations. | |
void | AddCrossSectionData (G4String materialName, G4String particleName, G4String fileCS, G4double scaleFactor) |
AddCrossSectionData Method used during the initialization of the model class to add a new material. It adds a material to the model and fills vectors with informations. Not every model needs differential cross sections. | |
void | LoadCrossSectionData (const G4String &particleName) |
LoadCrossSectionData Method to loop on all the registered materials in the model and load the corresponding data. | |
virtual void | ReadDiffCSFile (const G4String &materialName, const G4String &particleName, const G4String &path, const G4double scaleFactor) |
ReadDiffCSFile Virtual method that need to be implemented if one wish to use the differential cross sections. The read method for that kind of information is not standardized yet. | |
void | EnableForMaterialAndParticle (const G4String &materialName, const G4String &particleName) |
EnableMaterialAndParticle. | |
The G4VDNAModel class.
All the models using the DNA material management should inherit from that class. The goal is to allow the use of the material management system with little code interferences within the model classes.
Definition at line 49 of file G4VDNAModel.hh.
|
protected |
Definition at line 185 of file G4VDNAModel.hh.
|
protected |
Definition at line 184 of file G4VDNAModel.hh.
|
protected |
Definition at line 183 of file G4VDNAModel.hh.
G4VDNAModel Constructeur of the G4VDNAModel class.
nam | |
applyToMaterial |
Definition at line 35 of file G4VDNAModel.cc.
|
virtual |
~G4VDNAModel
Definition at line 41 of file G4VDNAModel.cc.
|
protected |
AddCrossSectionData Method used during the initialization of the model class to add a new material. It adds a material to the model and fills vectors with informations. Not every model needs differential cross sections.
materialName | |
particleName | |
fileCS | |
scaleFactor |
Definition at line 67 of file G4VDNAModel.cc.
|
protected |
AddCrossSectionData Method used during the initialization of the model class to add a new material. It adds a material to the model and fills vectors with informations.
materialName | |
particleName | |
fileCS | |
fileDiffCS | |
scaleFactor |
Definition at line 58 of file G4VDNAModel.cc.
Referenced by G4DNAPTBElasticModel::Initialise(), G4DNAPTBExcitationModel::Initialise(), and G4DNAPTBIonisationModel::Initialise().
BuildApplyToMatVect Build the material name vector which is used to know the materials the user want to include in the model.
materials |
Definition at line 139 of file G4VDNAModel.cc.
Referenced by LoadCrossSectionData().
|
pure virtual |
CrossSectionPerVolume Every model must implement its own CrossSectionPerVolume method. It is used by the process to determine the step path and must return a cross section times a number of molecules per volume unit.
material | |
materialName | |
p | |
ekin | |
emin | |
emax |
Implemented in G4DNADummyModel, G4DNAPTBElasticModel, G4DNAPTBExcitationModel, G4DNAPTBIonisationModel, and G4DNAVacuumModel.
Referenced by G4DNAModelInterface::CrossSectionPerVolume().
|
protected |
EnableMaterialAndParticle.
materialName | |
particleName | Meant to fill fTableData with 0 for the specified material and particle, therefore allowing the ModelInterface class to proceed with the material and particle even if no data are registered here. The data should obviously be registered somewhere in the child class. This method is here to allow an easy use of the no-ModelInterface dna models within the ModelInterface system. |
Definition at line 134 of file G4VDNAModel.cc.
Referenced by G4DNAVacuumModel::Initialise(), and G4DNADummyModel::Initialise().
GetHighEnergyLimit.
material | |
particle |
Definition at line 153 of file G4VDNAModel.hh.
Referenced by G4DNAPTBElasticModel::CrossSectionPerVolume(), G4DNAPTBExcitationModel::CrossSectionPerVolume(), G4DNAPTBIonisationModel::CrossSectionPerVolume(), G4DNAPTBElasticModel::SampleSecondaries(), and G4DNAPTBIonisationModel::SampleSecondaries().
GetLowEnergyLimit.
material | |
particle |
Definition at line 161 of file G4VDNAModel.hh.
Referenced by G4DNAPTBElasticModel::CrossSectionPerVolume(), G4DNAPTBExcitationModel::CrossSectionPerVolume(), G4DNAPTBIonisationModel::CrossSectionPerVolume(), G4DNAPTBElasticModel::SampleSecondaries(), and G4DNAPTBIonisationModel::SampleSecondaries().
|
inline |
GetName.
Definition at line 145 of file G4VDNAModel.hh.
Referenced by IsMaterialDefine().
|
inlineprotected |
GetTableData.
Definition at line 193 of file G4VDNAModel.hh.
Referenced by G4DNAPTBElasticModel::CrossSectionPerVolume(), G4DNAPTBExcitationModel::CrossSectionPerVolume(), G4DNAPTBIonisationModel::CrossSectionPerVolume(), and RandomSelectShell().
|
pure virtual |
Initialise Each model must implement an Initialize method.
particle | |
cuts |
Implemented in G4DNAVacuumModel, G4DNAPTBElasticModel, G4DNADummyModel, G4DNAPTBExcitationModel, and G4DNAPTBIonisationModel.
IsMaterialDefine Check if the given material is defined in the simulation.
materialName |
Definition at line 237 of file G4VDNAModel.cc.
IsMaterialExistingInModel Check if the given material is defined in the current model class.
materialName |
Definition at line 257 of file G4VDNAModel.cc.
Referenced by IsParticleExistingInModelForMaterial().
G4bool G4VDNAModel::IsParticleExistingInModelForMaterial | ( | const G4String & | particleName, |
const G4String & | materialName | ||
) |
IsParticleExistingInModelForMaterial To check two things: 1- is the material existing in model ? 2- if yes, is the particle defined for that material ?
particleName | |
materialName |
Definition at line 271 of file G4VDNAModel.cc.
|
protected |
LoadCrossSectionData Method to loop on all the registered materials in the model and load the corresponding data.
Definition at line 75 of file G4VDNAModel.cc.
Referenced by G4DNAPTBElasticModel::Initialise(), G4DNAPTBExcitationModel::Initialise(), and G4DNAPTBIonisationModel::Initialise().
|
protected |
RandomSelectShell Method to randomely select a shell from the data table uploaded. The size of the table (number of columns) is used to determine the total number of possible shells.
k | |
particle | |
materialName |
Definition at line 182 of file G4VDNAModel.cc.
Referenced by G4DNAPTBExcitationModel::SampleSecondaries(), and G4DNAPTBIonisationModel::SampleSecondaries().
|
protected |
ReadAndSaveCSFile Read and save a "simple" cross section file : use of G4DNACrossSectionDataSet->loadData()
materialName | |
particleName | |
file | |
scaleFactor |
Definition at line 174 of file G4VDNAModel.cc.
Referenced by LoadCrossSectionData().
|
protectedvirtual |
ReadDiffCSFile Virtual method that need to be implemented if one wish to use the differential cross sections. The read method for that kind of information is not standardized yet.
materialName | |
particleName | |
path | |
scaleFactor |
Definition at line 126 of file G4VDNAModel.cc.
Referenced by LoadCrossSectionData().
|
pure virtual |
SampleSecondaries Each model must implement SampleSecondaries to decide if a particle will be created after the ModelInterface or if any charateristic of the incident particle will change.
materialName | |
particleChangeForGamma | |
tmin | |
tmax |
Implemented in G4DNADummyModel, G4DNAPTBElasticModel, G4DNAPTBExcitationModel, G4DNAPTBIonisationModel, and G4DNAVacuumModel.
Referenced by G4DNAModelInterface::SampleSecondaries().
|
inline |
SetHighEnergyLimit.
material | |
particle | |
lim |
Definition at line 169 of file G4VDNAModel.hh.
Referenced by G4DNAPTBElasticModel::Initialise(), G4DNADummyModel::Initialise(), G4DNAPTBExcitationModel::Initialise(), and G4DNAPTBIonisationModel::Initialise().
|
inline |
SetLowEnergyLimit.
material | |
particle | |
lim |
Definition at line 177 of file G4VDNAModel.hh.
Referenced by G4DNAPTBElasticModel::Initialise(), G4DNADummyModel::Initialise(), G4DNAPTBExcitationModel::Initialise(), and G4DNAPTBIonisationModel::Initialise().