Geant4-11
|
The G4DNAPTBElasticModel class This class implements the elastic model for the DNA materials and precursors. More...
#include <G4DNAPTBElasticModel.hh>
Public Member Functions | |
virtual G4double | CrossSectionPerVolume (const G4Material *material, const G4String &materialName, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) |
CrossSectionPerVolume This method is mandatory for any model class. It finds and return the cross section value for the current material, particle and energy values. The number of molecule per volume is not used here but in the G4DNAModelInterface class. More... | |
G4DNAPTBElasticModel (const G4String &applyToMaterial="all", const G4ParticleDefinition *p=0, const G4String &nam="DNAPTBElasticModel") | |
G4DNAPTBElasticModel Constructor. More... | |
G4double | GetHighELimit (const G4String &material, const G4String &particle) |
GetHighEnergyLimit. More... | |
G4double | GetLowELimit (const G4String &material, const G4String &particle) |
GetLowEnergyLimit. More... | |
G4String | GetName () |
GetName. More... | |
virtual void | Initialise (const G4ParticleDefinition *particle, const G4DataVector &, G4ParticleChangeForGamma *fpChangeForGamme=nullptr) |
Initialise Mandatory method for every model class. The material/particle for which the model can be used have to be added here through the AddCrossSectionData method. Then the LoadCrossSectionData method must be called to trigger the load process. Scale factors to be applied to the cross section can be defined here. More... | |
G4bool | IsMaterialDefine (const G4String &materialName) |
IsMaterialDefine Check if the given material is defined in the simulation. More... | |
G4bool | IsMaterialExistingInModel (const G4String &materialName) |
IsMaterialExistingInModel Check if the given material is defined in the current model class. More... | |
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 ? More... | |
virtual void | SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4String &materialName, const G4DynamicParticle *, G4ParticleChangeForGamma *particleChangeForGamma, G4double tmin, G4double tmax) |
SampleSecondaries Method called after CrossSectionPerVolume if the process is the one which is selected (according to the sampling on the calculated path length). Here, the characteristics of the incident and created (if any) particle(s) are set (energy, momentum ...). More... | |
void | SetHighELimit (const G4String &material, const G4String &particle, G4double lim) |
SetHighEnergyLimit. More... | |
void | SetLowELimit (const G4String &material, const G4String &particle, G4double lim) |
SetLowEnergyLimit. More... | |
virtual | ~G4DNAPTBElasticModel () |
~G4DNAPTBElasticModel Destructor More... | |
Protected Types | |
typedef std::map< G4String, G4double >::const_iterator | ItCompoMapData |
typedef std::map< G4String, std::map< G4String, G4double > > | RatioMapData |
typedef std::map< G4String, std::map< G4String, G4DNACrossSectionDataSet *, std::less< G4String > > > | TableMapData |
Protected Member Functions | |
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. More... | |
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. More... | |
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. More... | |
void | EnableForMaterialAndParticle (const G4String &materialName, const G4String &particleName) |
EnableMaterialAndParticle. More... | |
TableMapData * | GetTableData () |
GetTableData. More... | |
void | LoadCrossSectionData (const G4String &particleName) |
LoadCrossSectionData Method to loop on all the registered materials in the model and load the corresponding data. More... | |
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. More... | |
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() More... | |
Private Types | |
typedef std::map< G4String, std::map< G4String, std::map< double, std::map< double, double > > > > | TriDimensionMap |
typedef std::map< G4String, std::map< G4String, std::map< double, std::vector< double > > > > | VecMap |
Private Member Functions | |
G4DNAPTBElasticModel (G4DNAPTBElasticModel &) | |
G4double | LinLinInterpolate (G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2) |
LinLinInterpolate. More... | |
G4double | LinLogInterpolate (G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2) |
LinLogInterpolate. More... | |
G4double | LogLogInterpolate (G4double e1, G4double e2, G4double e, G4double xs1, G4double xs2) |
LogLogInterpolate. More... | |
G4DNAPTBElasticModel & | operator= (const G4DNAPTBElasticModel &right) |
G4double | QuadInterpolator (G4double e11, G4double e12, G4double e21, G4double e22, G4double x11, G4double x12, G4double x21, G4double x22, G4double t1, G4double t2, G4double t, G4double e) |
QuadInterpolator. More... | |
G4double | RandomizeCosTheta (G4double k, const G4String &materialName) |
RandomizeCosTheta. More... | |
void | ReadDiffCSFile (const G4String &materialName, const G4String &particleName, const G4String &file, const G4double) |
ReadDiffCSFile Method to read the differential cross section files. This method is not standard yet so every model must implement its own. More... | |
G4double | Theta (G4ParticleDefinition *fParticleDefinition, G4double k, G4double integrDiff, const G4String &materialName) |
Theta To return an angular theta value from the differential file. This method uses interpolations to calculate the theta value. More... | |
Private Attributes | |
TriDimensionMap | diffCrossSectionData |
A map: [materialName][particleName]=DiffCrossSectionTable. More... | |
VecMap | eValuesVect |
std::map< G4String, std::map< G4String, G4double > > | fHighEnergyLimits |
List the high energy limits. More... | |
G4double | fKillBelowEnergy |
energy kill limit More... | |
std::map< G4String, std::map< G4String, G4double > > | fLowEnergyLimits |
List the low energy limits. More... | |
std::vector< G4String > | fModelCSFiles |
List the cross section data files. More... | |
std::vector< G4String > | fModelDiffCSFiles |
List the differential corss section data files. More... | |
std::vector< G4String > | fModelMaterials |
List the materials that can be activated (and will be by default) within the model. More... | |
std::vector< G4String > | fModelParticles |
List the particles that can be activated within the model. More... | |
std::vector< G4double > | fModelScaleFactors |
List the model scale factors (they could change with material) More... | |
G4String | fName |
model name More... | |
const G4String | fStringOfMaterials |
fStringOfMaterials The user can decide to specify by hand which are the materials the be activated among those implemented in the model. If the user does then only the specified materials contained in this string variable will be activated. The string is like: mat1/mat2/mat3/mat4 More... | |
TableMapData | fTableData |
fTableData It contains the cross section data and can be used like: dataTable=fTableData[material][particle] More... | |
std::map< G4String, double > | killBelowEnergyTable |
map to save the different energy kill limits for the materials More... | |
std::map< G4String, std::map< G4String, std::vector< double > > > | tValuesVec |
map with vectors containing all the incident (T) energy of the differential file More... | |
G4int | verboseLevel |
verbose level More... | |
The G4DNAPTBElasticModel class This class implements the elastic model for the DNA materials and precursors.
Definition at line 48 of file G4DNAPTBElasticModel.hh.
|
protectedinherited |
Definition at line 185 of file G4VDNAModel.hh.
|
protectedinherited |
Definition at line 184 of file G4VDNAModel.hh.
|
protectedinherited |
Definition at line 183 of file G4VDNAModel.hh.
|
private |
Definition at line 125 of file G4DNAPTBElasticModel.hh.
|
private |
Definition at line 128 of file G4DNAPTBElasticModel.hh.
G4DNAPTBElasticModel::G4DNAPTBElasticModel | ( | const G4String & | applyToMaterial = "all" , |
const G4ParticleDefinition * | p = 0 , |
||
const G4String & | nam = "DNAPTBElasticModel" |
||
) |
G4DNAPTBElasticModel Constructor.
applyToMaterial | |
p | |
nam |
Definition at line 38 of file G4DNAPTBElasticModel.cc.
References eV, fKillBelowEnergy, G4cout, G4endl, and verboseLevel.
|
virtual |
|
private |
|
protectedinherited |
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.
References G4VDNAModel::fModelCSFiles, G4VDNAModel::fModelMaterials, G4VDNAModel::fModelParticles, and G4VDNAModel::fModelScaleFactors.
|
protectedinherited |
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.
References G4VDNAModel::fModelCSFiles, G4VDNAModel::fModelDiffCSFiles, G4VDNAModel::fModelMaterials, G4VDNAModel::fModelParticles, and G4VDNAModel::fModelScaleFactors.
Referenced by Initialise(), G4DNAPTBExcitationModel::Initialise(), and G4DNAPTBIonisationModel::Initialise().
|
protectedinherited |
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 G4VDNAModel::LoadCrossSectionData().
|
virtual |
CrossSectionPerVolume This method is mandatory for any model class. It finds and return the cross section value for the current material, particle and energy values. The number of molecule per volume is not used here but in the G4DNAModelInterface class.
material | |
materialName | |
p | |
ekin | |
emin | |
emax |
Implements G4VDNAModel.
Definition at line 295 of file G4DNAPTBElasticModel.cc.
References cm, DBL_MAX, eV, fKillBelowEnergy, G4cout, G4endl, G4VDNAModel::GetHighELimit(), G4VDNAModel::GetLowELimit(), G4ParticleDefinition::GetParticleName(), G4VDNAModel::GetTableData(), and verboseLevel.
|
protectedinherited |
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.
References G4VDNAModel::fTableData.
Referenced by G4DNAVacuumModel::Initialise(), and G4DNADummyModel::Initialise().
|
inlineinherited |
GetHighEnergyLimit.
material | |
particle |
Definition at line 153 of file G4VDNAModel.hh.
References G4VDNAModel::fHighEnergyLimits, and eplot::material.
Referenced by CrossSectionPerVolume(), G4DNAPTBExcitationModel::CrossSectionPerVolume(), G4DNAPTBIonisationModel::CrossSectionPerVolume(), SampleSecondaries(), and G4DNAPTBIonisationModel::SampleSecondaries().
|
inlineinherited |
GetLowEnergyLimit.
material | |
particle |
Definition at line 161 of file G4VDNAModel.hh.
References G4VDNAModel::fLowEnergyLimits, and eplot::material.
Referenced by CrossSectionPerVolume(), G4DNAPTBExcitationModel::CrossSectionPerVolume(), G4DNAPTBIonisationModel::CrossSectionPerVolume(), SampleSecondaries(), and G4DNAPTBIonisationModel::SampleSecondaries().
|
inlineinherited |
GetName.
Definition at line 145 of file G4VDNAModel.hh.
References G4VDNAModel::fName.
Referenced by G4VDNAModel::IsMaterialDefine().
|
inlineprotectedinherited |
GetTableData.
Definition at line 193 of file G4VDNAModel.hh.
References G4VDNAModel::fTableData.
Referenced by CrossSectionPerVolume(), G4DNAPTBExcitationModel::CrossSectionPerVolume(), G4DNAPTBIonisationModel::CrossSectionPerVolume(), and G4VDNAModel::RandomSelectShell().
|
virtual |
Initialise Mandatory method for every model class. The material/particle for which the model can be used have to be added here through the AddCrossSectionData method. Then the LoadCrossSectionData method must be called to trigger the load process. Scale factors to be applied to the cross section can be defined here.
Implements G4VDNAModel.
Definition at line 67 of file G4DNAPTBElasticModel.cc.
References G4VDNAModel::AddCrossSectionData(), cm, G4Electron::ElectronDefinition(), eV, G4cout, G4endl, G4ParticleDefinition::GetParticleName(), keV, G4VDNAModel::LoadCrossSectionData(), G4VDNAModel::SetHighELimit(), G4VDNAModel::SetLowELimit(), and verboseLevel.
IsMaterialDefine Check if the given material is defined in the simulation.
materialName |
Definition at line 237 of file G4VDNAModel.cc.
References G4Material::GetMaterialTable(), and G4VDNAModel::GetName().
IsMaterialExistingInModel Check if the given material is defined in the current model class.
materialName |
Definition at line 257 of file G4VDNAModel.cc.
References G4VDNAModel::fTableData.
Referenced by G4VDNAModel::IsParticleExistingInModelForMaterial().
|
inherited |
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.
References G4VDNAModel::fTableData, and G4VDNAModel::IsMaterialExistingInModel().
|
private |
LinLinInterpolate.
e1 | |
e2 | |
e | |
xs1 | |
xs2 |
Definition at line 470 of file G4DNAPTBElasticModel.cc.
References d1, d2, e1, and e2.
Referenced by QuadInterpolator().
|
private |
LinLogInterpolate.
e1 | |
e2 | |
e | |
xs1 | |
xs2 |
Definition at line 456 of file G4DNAPTBElasticModel.cc.
|
protectedinherited |
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.
References G4VDNAModel::BuildApplyToMatVect(), FatalException, G4VDNAModel::fModelCSFiles, G4VDNAModel::fModelDiffCSFiles, G4VDNAModel::fModelMaterials, G4VDNAModel::fModelParticles, G4VDNAModel::fModelScaleFactors, G4VDNAModel::fStringOfMaterials, G4Exception(), G4VDNAModel::ReadAndSaveCSFile(), and G4VDNAModel::ReadDiffCSFile().
Referenced by Initialise(), G4DNAPTBExcitationModel::Initialise(), and G4DNAPTBIonisationModel::Initialise().
|
private |
LogLogInterpolate.
e1 | |
e2 | |
e | |
xs1 | |
xs2 |
Definition at line 484 of file G4DNAPTBElasticModel.cc.
|
private |
|
private |
QuadInterpolator.
e11 | |
e12 | |
e21 | |
e22 | |
x11 | |
x12 | |
x21 | |
x22 | |
t1 | |
t2 | |
t | |
e |
Definition at line 499 of file G4DNAPTBElasticModel.cc.
References LinLinInterpolate().
Referenced by Theta().
|
private |
RandomizeCosTheta.
k | |
materialName |
Definition at line 529 of file G4DNAPTBElasticModel.cc.
References G4Electron::ElectronDefinition(), eV, G4UniformRand, pi, Theta(), and anonymous_namespace{G4HadPhaseSpaceGenbod.cc}::uniformRand().
Referenced by SampleSecondaries().
|
protectedinherited |
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.
References FatalException, G4VEMDataSet::FindValue(), G4Exception(), G4UniformRand, G4DNACrossSectionDataSet::GetComponent(), G4VDNAModel::GetTableData(), CLHEP::detail::n, G4DNACrossSectionDataSet::NumberOfComponents(), and pos.
Referenced by G4DNAPTBExcitationModel::SampleSecondaries(), and G4DNAPTBIonisationModel::SampleSecondaries().
|
protectedinherited |
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.
References eV, geant4_check_module_cycles::file, and G4VDNAModel::fTableData.
Referenced by G4VDNAModel::LoadCrossSectionData().
|
privatevirtual |
ReadDiffCSFile Method to read the differential cross section files. This method is not standard yet so every model must implement its own.
materialName | |
particleName | |
file |
Reimplemented from G4VDNAModel.
Definition at line 197 of file G4DNAPTBElasticModel.cc.
References diffCrossSectionData, eValuesVect, FatalException, geant4_check_module_cycles::file, G4Exception(), mcscore::test(), and tValuesVec.
|
virtual |
SampleSecondaries Method called after CrossSectionPerVolume if the process is the one which is selected (according to the sampling on the calculated path length). Here, the characteristics of the incident and created (if any) particle(s) are set (energy, momentum ...).
materialName | |
particleChangeForGamma | |
tmin | |
tmax |
Implements G4VDNAModel.
Definition at line 346 of file G4DNAPTBElasticModel.cc.
References CLHEP::Hep3Vector::cross(), fKillBelowEnergy, fStopAndKill, G4cout, G4endl, G4UniformRand, G4VDNAModel::GetHighELimit(), G4DynamicParticle::GetKineticEnergy(), G4VDNAModel::GetLowELimit(), G4DynamicParticle::GetMomentumDirection(), G4DynamicParticle::GetParticleDefinition(), G4ParticleDefinition::GetParticleName(), CLHEP::Hep3Vector::orthogonal(), pi, G4VParticleChange::ProposeLocalEnergyDeposit(), G4ParticleChangeForGamma::ProposeMomentumDirection(), G4VParticleChange::ProposeTrackStatus(), RandomizeCosTheta(), G4ParticleChangeForGamma::SetProposedKineticEnergy(), CLHEP::Hep3Vector::unit(), and verboseLevel.
|
inlineinherited |
SetHighEnergyLimit.
material | |
particle | |
lim |
Definition at line 169 of file G4VDNAModel.hh.
References G4VDNAModel::fHighEnergyLimits, and eplot::material.
Referenced by Initialise(), G4DNADummyModel::Initialise(), G4DNAPTBExcitationModel::Initialise(), and G4DNAPTBIonisationModel::Initialise().
|
inlineinherited |
SetLowEnergyLimit.
material | |
particle | |
lim |
Definition at line 177 of file G4VDNAModel.hh.
References G4VDNAModel::fLowEnergyLimits, and eplot::material.
Referenced by Initialise(), G4DNADummyModel::Initialise(), G4DNAPTBExcitationModel::Initialise(), and G4DNAPTBIonisationModel::Initialise().
|
private |
Theta To return an angular theta value from the differential file. This method uses interpolations to calculate the theta value.
fParticleDefinition | |
k | |
integrDiff | |
materialName |
Definition at line 402 of file G4DNAPTBElasticModel.cc.
References diffCrossSectionData, G4Electron::ElectronDefinition(), eValuesVect, G4ParticleDefinition::GetParticleName(), QuadInterpolator(), and tValuesVec.
Referenced by RandomizeCosTheta().
|
private |
A map: [materialName][particleName]=DiffCrossSectionTable.
Definition at line 126 of file G4DNAPTBElasticModel.hh.
Referenced by ReadDiffCSFile(), and Theta().
|
private |
map with vectors containing all the output energy (E) of the differential file
Definition at line 129 of file G4DNAPTBElasticModel.hh.
Referenced by ReadDiffCSFile(), and Theta().
List the high energy limits.
Definition at line 301 of file G4VDNAModel.hh.
Referenced by G4VDNAModel::GetHighELimit(), and G4VDNAModel::SetHighELimit().
|
private |
energy kill limit
Definition at line 123 of file G4DNAPTBElasticModel.hh.
Referenced by CrossSectionPerVolume(), G4DNAPTBElasticModel(), and SampleSecondaries().
List the low energy limits.
Definition at line 300 of file G4VDNAModel.hh.
Referenced by G4VDNAModel::GetLowELimit(), and G4VDNAModel::SetLowELimit().
|
privateinherited |
List the cross section data files.
Definition at line 296 of file G4VDNAModel.hh.
Referenced by G4VDNAModel::AddCrossSectionData(), and G4VDNAModel::LoadCrossSectionData().
|
privateinherited |
List the differential corss section data files.
Definition at line 297 of file G4VDNAModel.hh.
Referenced by G4VDNAModel::AddCrossSectionData(), and G4VDNAModel::LoadCrossSectionData().
|
privateinherited |
List the materials that can be activated (and will be by default) within the model.
Definition at line 294 of file G4VDNAModel.hh.
Referenced by G4VDNAModel::AddCrossSectionData(), and G4VDNAModel::LoadCrossSectionData().
|
privateinherited |
List the particles that can be activated within the model.
Definition at line 295 of file G4VDNAModel.hh.
Referenced by G4VDNAModel::AddCrossSectionData(), and G4VDNAModel::LoadCrossSectionData().
|
privateinherited |
List the model scale factors (they could change with material)
Definition at line 298 of file G4VDNAModel.hh.
Referenced by G4VDNAModel::AddCrossSectionData(), and G4VDNAModel::LoadCrossSectionData().
|
privateinherited |
|
privateinherited |
fStringOfMaterials The user can decide to specify by hand which are the materials the be activated among those implemented in the model. If the user does then only the specified materials contained in this string variable will be activated. The string is like: mat1/mat2/mat3/mat4
Definition at line 286 of file G4VDNAModel.hh.
Referenced by G4VDNAModel::LoadCrossSectionData().
|
privateinherited |
fTableData It contains the cross section data and can be used like: dataTable=fTableData[material][particle]
Definition at line 292 of file G4VDNAModel.hh.
Referenced by G4VDNAModel::EnableForMaterialAndParticle(), G4VDNAModel::GetTableData(), G4VDNAModel::IsMaterialExistingInModel(), G4VDNAModel::IsParticleExistingInModelForMaterial(), G4VDNAModel::ReadAndSaveCSFile(), and G4VDNAModel::~G4VDNAModel().
|
private |
map to save the different energy kill limits for the materials
Definition at line 122 of file G4DNAPTBElasticModel.hh.
|
private |
map with vectors containing all the incident (T) energy of the differential file
Definition at line 130 of file G4DNAPTBElasticModel.hh.
Referenced by ReadDiffCSFile(), and Theta().
|
private |
verbose level
Definition at line 121 of file G4DNAPTBElasticModel.hh.
Referenced by CrossSectionPerVolume(), G4DNAPTBElasticModel(), Initialise(), and SampleSecondaries().