Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4DNADummyModel Class Reference

#include <G4DNADummyModel.hh>

+ Inheritance diagram for G4DNADummyModel:

Public Member Functions

 G4DNADummyModel (const G4String &applyToMaterial, const G4ParticleDefinition *p, const G4String &nam, G4VEmModel *emModel)
 
 ~G4DNADummyModel ()
 
virtual void Initialise (const G4ParticleDefinition *particle, const G4DataVector &= *(new G4DataVector()), G4ParticleChangeForGamma *changeForGamme=nullptr)
 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)
 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, G4double tmax)
 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.
 
const G4VEmModelGetEmModel () const
 
G4VEmModelGetEmModel ()
 
- Public Member Functions inherited from G4VDNAModel
 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.
 

Additional Inherited Members

- Protected Types inherited from G4VDNAModel
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 inherited from G4VDNAModel
TableMapDataGetTableData ()
 GetTableData.
 
std::vector< G4StringBuildApplyToMatVect (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.
 

Detailed Description

Definition at line 40 of file G4DNADummyModel.hh.

Constructor & Destructor Documentation

◆ G4DNADummyModel()

G4DNADummyModel::G4DNADummyModel ( const G4String applyToMaterial,
const G4ParticleDefinition p,
const G4String nam,
G4VEmModel emModel 
)

Definition at line 35 of file G4DNADummyModel.cc.

36 : G4VDNAModel(nam, applyToMaterial)
37{
38 fpEmModel = emModel;
39 fpParticleDef = p;
40}
The G4VDNAModel class.
Definition: G4VDNAModel.hh:50

◆ ~G4DNADummyModel()

G4DNADummyModel::~G4DNADummyModel ( )

Definition at line 42 of file G4DNADummyModel.cc.

43{
44 // There is no need to delete the model because it will be done in some G4 class.
45 //if(fpEmModel) delete fpEmModel;
46}

Member Function Documentation

◆ CrossSectionPerVolume()

G4double G4DNADummyModel::CrossSectionPerVolume ( const G4Material material,
const G4String materialName,
const G4ParticleDefinition p,
G4double  ekin,
G4double  emin,
G4double  emax 
)
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.

Parameters
material
materialName
p
ekin
emin
emax
Returns
crossSection*numberOfMoleculesPerVolumeUnit

Implements G4VDNAModel.

Definition at line 61 of file G4DNADummyModel.cc.

62{
63 G4double crossSectionTimesDensity = fpEmModel->CrossSectionPerVolume(material, p, ekin, emin, emax);
64 G4double crossSection = crossSectionTimesDensity / GetNumMoleculePerVolumeUnitForMaterial(G4Material::GetMaterial("G4_WATER") );
65
66 return crossSection;
67}
double G4double
Definition: G4Types.hh:83
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:691
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:181

◆ GetEmModel() [1/2]

G4VEmModel * G4DNADummyModel::GetEmModel ( )
inline

Definition at line 67 of file G4DNADummyModel.hh.

67{return fpEmModel;}

◆ GetEmModel() [2/2]

const G4VEmModel * G4DNADummyModel::GetEmModel ( ) const
inline

Definition at line 66 of file G4DNADummyModel.hh.

66{return fpEmModel;}

◆ Initialise()

void G4DNADummyModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector cuts = *(new G4DataVector()),
G4ParticleChangeForGamma fpChangeForGamme = nullptr 
)
virtual

Initialise Each model must implement an Initialize method.

Parameters
particle
cuts

Implements G4VDNAModel.

Definition at line 48 of file G4DNADummyModel.cc.

49{
51
52 fpEmModel->SetParticleChange(changeForGamme, nullptr);
53 fpEmModel->Initialise(particle, v);
54
55 // MatManagSys
56 EnableForMaterialAndParticle("G4_WATER", fpParticleDef->GetParticleName() );
57 SetLowELimit("G4_WATER", fpParticleDef->GetParticleName(), fpEmModel->LowEnergyLimit() );
58 SetHighELimit("G4_WATER",fpParticleDef->GetParticleName(), fpEmModel->HighEnergyLimit() );
59}
const std::vector< G4double > * GetNumMolPerVolTableFor(const G4Material *) const
Retrieve a table of molecular densities (number of molecules per unit volume) in the G4 unit system f...
static G4DNAMolecularMaterial * Instance()
const G4String & GetParticleName() const
void SetHighELimit(const G4String &material, const G4String &particle, G4double lim)
SetHighEnergyLimit.
Definition: G4VDNAModel.hh:169
void SetLowELimit(const G4String &material, const G4String &particle, G4double lim)
SetLowEnergyLimit.
Definition: G4VDNAModel.hh:177
void EnableForMaterialAndParticle(const G4String &materialName, const G4String &particleName)
EnableMaterialAndParticle.
Definition: G4VDNAModel.cc:134
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
void SetParticleChange(G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
Definition: G4VEmModel.cc:390
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)=0

◆ SampleSecondaries()

void G4DNADummyModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  ,
const G4MaterialCutsCouple ,
const G4String materialName,
const G4DynamicParticle ,
G4ParticleChangeForGamma particleChangeForGamma,
G4double  tmin,
G4double  tmax 
)
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.

Parameters
materialName
particleChangeForGamma
tmin
tmax

Implements G4VDNAModel.

Definition at line 69 of file G4DNADummyModel.cc.

70{
71 fpEmModel->SampleSecondaries(a, b, c, tmin, tmax);
72}
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0

The documentation for this class was generated from the following files: