Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4VDNAModel Class Referenceabstract

The G4VDNAModel class. More...

#include <G4VDNAModel.hh>

+ Inheritance diagram for G4VDNAModel:

Public Member Functions

 G4VDNAModel (const G4String &nam, const G4String &applyToMaterial="")
 G4VDNAModel Constructeur of the G4VDNAModel class.
 
 ~G4VDNAModel () override
 ~G4VDNAModel
 
void Initialise (const G4ParticleDefinition *particle, const G4DataVector &cuts) override=0
 Initialise Each model must implement an Initialize method.
 
G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override=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.
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0, G4double tmax=DBL_MAX) override=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 size_t &materialID)
 IsMaterialDefine Check if the given material is defined in the simulation.
 
G4bool IsParticleExistingInModelForMaterial (const G4ParticleDefinition *particleName, const size_t &materialID)
 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 size_t &materialID, const G4ParticleDefinition *particle)
 GetHighEnergyLimit.
 
G4double GetLowELimit (const size_t &materialID, const G4ParticleDefinition *particle)
 GetLowEnergyLimit.
 
void SetHighELimit (const size_t &materialID, const G4ParticleDefinition *particle, G4double lim)
 SetHighEnergyLimit.
 
void SetLowELimit (const size_t &materialID, const G4ParticleDefinition *particle, G4double lim)
 SetLowEnergyLimit.
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, const G4double &length, G4double &eloss)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
virtual void FillNumberOfSecondaries (G4int &numberOfTriplets, G4int &numberOfRecoil)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementGetCurrentElement (const G4Material *mat=nullptr) const
 
G4int SelectRandomAtomNumber (const G4Material *) const
 
const G4IsotopeGetCurrentIsotope (const G4Element *elm=nullptr) const
 
G4int SelectIsotopeNumber (const G4Element *) const
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
void SetLPMFlag (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Types

using MaterialParticleMapData
 

Protected Member Functions

MaterialParticleMapDataGetData ()
 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 size_t &materialID, const G4ParticleDefinition *p, const G4String &file, const G4double &scaleFactor)
 ReadAndSaveCSFile Read and save a "simple" cross section file : use of G4DNACrossSectionDataSet->loadData()
 
G4int RandomSelectShell (const G4double &k, const G4ParticleDefinition *particle, const size_t &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 (const size_t &materialName, const G4ParticleDefinition *particleName, const G4String &fileCS, const G4String &fileDiffCS, const 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 (const size_t &materialName, const G4ParticleDefinition *particleName, const G4String &fileCS, const 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 G4ParticleDefinition *particleName)
 LoadCrossSectionData Method to loop on all the registered materials in the model and load the corresponding data.
 
virtual void ReadDiffCSFile (const size_t &materialName, const G4ParticleDefinition *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 size_t &materialID, const G4ParticleDefinition *p)
 EnableMaterialAndParticle.
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Protected Attributes

G4bool isInitialised = false
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4PhysicsTablexSectionTable = nullptr
 
const G4MaterialpBaseMaterial = nullptr
 
const std::vector< G4double > * theDensityFactor = nullptr
 
const std::vector< G4int > * theDensityIdx = nullptr
 
G4double inveplus
 
G4double pFactor = 1.0
 
std::size_t currentCoupleIndex = 0
 
std::size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

Detailed Description

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 50 of file G4VDNAModel.hh.

Member Typedef Documentation

◆ MaterialParticleMapData

Initial value:
std::map<size_t ,
std::map<const G4ParticleDefinition*, std::unique_ptr<G4DNACrossSectionDataSet>>>

Definition at line 173 of file G4VDNAModel.hh.

Constructor & Destructor Documentation

◆ G4VDNAModel()

G4VDNAModel::G4VDNAModel ( const G4String & nam,
const G4String & applyToMaterial = "" )

G4VDNAModel Constructeur of the G4VDNAModel class.

Parameters
nam
applyToMaterial

Definition at line 36 of file G4VDNAModel.cc.

37 : G4VEmModel(nam), fStringOfMaterials(applyToMaterial)
38{}
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67

◆ ~G4VDNAModel()

G4VDNAModel::~G4VDNAModel ( )
overridedefault

~G4VDNAModel

Member Function Documentation

◆ AddCrossSectionData() [1/2]

void G4VDNAModel::AddCrossSectionData ( const size_t & materialName,
const G4ParticleDefinition * particleName,
const G4String & fileCS,
const G4double & scaleFactor )
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.

Parameters
materialName

param particleName

Parameters
fileCS

param scaleFactor

◆ AddCrossSectionData() [2/2]

void G4VDNAModel::AddCrossSectionData ( const size_t & materialName,
const G4ParticleDefinition * particleName,
const G4String & fileCS,
const G4String & fileDiffCS,
const G4double & scaleFactor )
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.

Parameters
materialName

param particleName

Parameters
fileCS

param fileDiffCS

Parameters
scaleFactor

Referenced by G4DNACPA100ElasticModel::Initialise(), G4DNACPA100ExcitationModel::Initialise(), G4DNACPA100IonisationModel::Initialise(), G4DNAPTBElasticModel::Initialise(), G4DNAPTBExcitationModel::Initialise(), and G4DNAPTBIonisationModel::Initialise().

◆ BuildApplyToMatVect()

std::vector< G4String > G4VDNAModel::BuildApplyToMatVect ( const G4String & materials)
protected

BuildApplyToMatVect Build the material name vector which is used to know the materials the user want to include in the model.

Parameters
materials
Returns
a vector with all the material names

Definition at line 138 of file G4VDNAModel.cc.

139{
140 // output material vector
141 std::vector<G4String> materialVect;
142
143 // if we don't find any "/" then it means we only have one "material" (could be the "all" option)
144 if (materials.find('/') == std::string::npos) {
145 // we add the material to the output vector
146 materialVect.push_back(materials);
147 }
148 // if we have several materials listed in the string then we must retrieve them
149 else {
150 G4String materialsNonIdentified = materials;
151
152 while (materialsNonIdentified.find_first_of('/') != std::string::npos) {
153 // we select the first material and stop at the "/" caracter
154 G4String mat = materialsNonIdentified.substr(0, materialsNonIdentified.find_first_of('/'));
155 materialVect.push_back(mat);
156
157 // we remove the previous material from the materialsNonIdentified string
158 materialsNonIdentified = materialsNonIdentified.substr(
159 materialsNonIdentified.find_first_of('/') + 1,
160 materialsNonIdentified.size() - materialsNonIdentified.find_first_of('/'));
161 }
162
163 // we don't find "/" anymore, it means we only have one material string left
164 // we get it
165 materialVect.push_back(materialsNonIdentified);
166 }
167
168 return materialVect;
169}

Referenced by LoadCrossSectionData().

◆ CrossSectionPerVolume()

G4double G4VDNAModel::CrossSectionPerVolume ( const G4Material * material,
const G4ParticleDefinition * p,
G4double ekin,
G4double emin,
G4double emax )
overridepure 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

param materialName

Parameters
p

param ekin

Parameters
emin
emax
Returns
crossSection*numberOfMoleculesPerVolumeUnit

Reimplemented from G4VEmModel.

Implemented in G4DNACPA100ElasticModel, G4DNACPA100ExcitationModel, G4DNACPA100IonisationModel, G4DNAPTBElasticModel, G4DNAPTBExcitationModel, G4DNAPTBIonisationModel, and G4DNAVacuumModel.

◆ EnableForMaterialAndParticle()

void G4VDNAModel::EnableForMaterialAndParticle ( const size_t & materialID,
const G4ParticleDefinition * p )
protected

EnableMaterialAndParticle.

Parameters
materialName
particleNameMeant 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 132 of file G4VDNAModel.cc.

134{
135 fData[materialID][p] = nullptr;
136}

Referenced by G4DNAVacuumModel::Initialise().

◆ GetData()

MaterialParticleMapData * G4VDNAModel::GetData ( )
inlineprotected

GetTableData.

Returns
a pointer to a map with the following structure: [materialName][particleName]=G4DNACrossSectionDataSet*

Definition at line 183 of file G4VDNAModel.hh.

184 {
185 return &fData;
186 }

Referenced by G4DNACPA100ElasticModel::CrossSectionPerVolume(), G4DNACPA100ExcitationModel::CrossSectionPerVolume(), G4DNACPA100IonisationModel::CrossSectionPerVolume(), G4DNAPTBElasticModel::CrossSectionPerVolume(), G4DNAPTBExcitationModel::CrossSectionPerVolume(), and G4DNAPTBIonisationModel::CrossSectionPerVolume().

◆ GetHighELimit()

◆ GetLowELimit()

◆ GetName()

◆ Initialise()

void G4VDNAModel::Initialise ( const G4ParticleDefinition * particle,
const G4DataVector & cuts )
overridepure virtual

Initialise Each model must implement an Initialize method.

Parameters
particle
cuts

Implements G4VEmModel.

Implemented in G4DNACPA100ElasticModel, G4DNACPA100ExcitationModel, G4DNACPA100IonisationModel, G4DNAPTBElasticModel, G4DNAPTBExcitationModel, G4DNAPTBIonisationModel, and G4DNAVacuumModel.

◆ IsMaterialDefine()

G4bool G4VDNAModel::IsMaterialDefine ( const size_t & materialID)

IsMaterialDefine Check if the given material is defined in the simulation.

Parameters
materialName
Returns
true if the material is defined in the simulation

Definition at line 228 of file G4VDNAModel.cc.

229{
230 // Check if the given material is defined in the simulation
231
232 G4bool exist(false);
233
234 G4double matTableSize = G4Material::GetMaterialTable()->size();
235
236 for (int i = 0; i < matTableSize; i++) {
237 if (materialID == G4Material::GetMaterialTable()->at(i)->GetIndex()) {
238 exist = true;
239 return exist;
240 }
241 }
242
243 G4Exception("G4VDNAModel::IsMaterialDefine", "em0003", FatalException,
244 "Materials are not defined!!");
245 return exist;
246}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
static G4MaterialTable * GetMaterialTable()

◆ IsParticleExistingInModelForMaterial()

G4bool G4VDNAModel::IsParticleExistingInModelForMaterial ( const G4ParticleDefinition * particleName,
const size_t & materialID )

IsParticleExistingInModelForMaterial To check two things: 1- is the material existing in model ? 2- if yes, is the particle defined for that material ?

Parameters
particleName
materialName
Returns
true if the particle/material couple is defined in the model

Definition at line 260 of file G4VDNAModel.cc.

262{
263 // To check two things:
264 // 1- is the material existing in model ?
265 // 2- if yes, is the particle defined for that material ?
266
267 if (IsMaterialExistingInModel(materialID)) {
268 for (const auto& it : fModelParticles) {
269 if (it == particleName) {
270 return true;
271 }
272 }
273 }
274 return false;
275}

◆ LoadCrossSectionData()

void G4VDNAModel::LoadCrossSectionData ( const G4ParticleDefinition * particleName)
protected

LoadCrossSectionData Method to loop on all the registered materials in the model and load the corresponding data.

Definition at line 64 of file G4VDNAModel.cc.

65{
66 G4String fileElectron, fileDiffElectron = "";
67 G4String materialName, modelParticleName;
68 G4double scaleFactor;
69 std::size_t materialID;
70
71 const G4ParticleDefinition* pParticle;
72
73 // construct applyToMatVect with materials specified by the user
74 std::vector<G4String> applyToMatVect = BuildApplyToMatVect(fStringOfMaterials);
75
76 // iterate on each material contained into the fStringOfMaterials variable (through
77 // applyToMatVect)
78 for (const auto & i : applyToMatVect) {
79 auto pMat = G4Material::GetMaterial(i, false);
80 if (i != "all" && pMat == nullptr) {
81 continue;
82 }
83
84 // We have selected a material coming from applyToMatVect
85 // We try to find if this material correspond to a model registered material
86 // If it is, then isMatFound becomes true
87 G4bool isMatFound = false;
88
89 // We iterate on each model registered materials to load the CS data
90 // We have to do a for loop because of the "all" option
91 // applyToMatVect[i] == "all" implies applyToMatVect.size()=1 and we want to iterate on all
92 // registered materials
93 for (std::size_t j = 0; j < fModelMaterials.size(); ++j) {
94 if (i == "all" || pMat->GetIndex() == fModelMaterials[j]) {
95 isMatFound = true;
96 materialID = fModelMaterials[j];
97 pParticle = fModelParticles[j];
98 fileElectron = fModelCSFiles[j];
99 if (!fModelDiffCSFiles.empty()) fileDiffElectron = fModelDiffCSFiles[j];
100 scaleFactor = fModelScaleFactors[j];
101
102 ReadAndSaveCSFile(materialID, pParticle, fileElectron, scaleFactor);
103
104 if (!fileDiffElectron.empty())
105 ReadDiffCSFile(materialID, pParticle, fileDiffElectron, scaleFactor);
106 }
107 }
108
109 // check if we found a correspondance, if not: fatal error
110 if (!isMatFound) {
111 std::ostringstream oss;
112 oss << i
113 << " material was not found. It means the material specified in the UserPhysicsList is "
114 "not a model material for ";
115 oss << particleName;
116 G4Exception("G4VDNAModel::LoadCrossSectionData", "em0003", FatalException, oss.str().c_str());
117 return;
118 }
119 }
120}
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
std::vector< G4String > BuildApplyToMatVect(const G4String &materials)
BuildApplyToMatVect Build the material name vector which is used to know the materials the user want ...
virtual void ReadDiffCSFile(const size_t &materialName, const G4ParticleDefinition *particleName, const G4String &path, const G4double &scaleFactor)
ReadDiffCSFile Virtual method that need to be implemented if one wish to use the differential cross s...
void ReadAndSaveCSFile(const size_t &materialID, const G4ParticleDefinition *p, const G4String &file, const G4double &scaleFactor)
ReadAndSaveCSFile Read and save a "simple" cross section file : use of G4DNACrossSectionDataSet->load...

Referenced by G4DNACPA100ElasticModel::Initialise(), G4DNACPA100ExcitationModel::Initialise(), G4DNACPA100IonisationModel::Initialise(), G4DNAPTBElasticModel::Initialise(), G4DNAPTBExcitationModel::Initialise(), and G4DNAPTBIonisationModel::Initialise().

◆ RandomSelectShell()

G4int G4VDNAModel::RandomSelectShell ( const G4double & k,
const G4ParticleDefinition * particle,
const size_t & materialName )
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.

Parameters
k

param particle

Parameters
materialName
Returns
the selected shell

Definition at line 179 of file G4VDNAModel.cc.

181{
182 G4int level = 0;
183
184 auto pos = fData[materialID].find(particle);
185
186 if (pos != fData[materialID].end()) {
187 G4DNACrossSectionDataSet* table = pos->second.get();
188
189 if (table != nullptr) {
190 auto valuesBuffer = new G4double[table->NumberOfComponents()];
191 auto n = (G4int)table->NumberOfComponents();
192 G4int i(n);
193 G4double value = 0.;
194
195 while (i > 0) {
196 --i;
197 valuesBuffer[i] = table->GetComponent(i)->FindValue(k);
198 value += valuesBuffer[i];
199 }
200
201 value *= G4UniformRand();
202
203 i = n;
204
205 while (i > 0) {
206 --i;
207
208 if (valuesBuffer[i] > value) {
209 delete[] valuesBuffer;
210 return i;
211 }
212 value -= valuesBuffer[i];
213 }
214
215 delete[] valuesBuffer;
216 }
217 }
218 else {
219 G4cout << "particle : " << particle->GetParticleName()
220 << " Materials : " << (*G4Material::GetMaterialTable())[materialID]->GetName() << " "
221 << this->GetName() << G4endl;
222 G4Exception("G4VDNAModel::RandomSelectShell", "em0002", FatalException,
223 "Model not applicable to particle type : ");
224 }
225 return level;
226}
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
#define G4UniformRand()
Definition Randomize.hh:52
size_t NumberOfComponents() const override
const G4VEMDataSet * GetComponent(G4int componentId) const override
const G4String & GetParticleName() const
G4String GetName()
GetName.
virtual G4double FindValue(G4double x, G4int componentId=0) const =0

Referenced by G4DNACPA100ExcitationModel::SampleSecondaries(), G4DNACPA100IonisationModel::SampleSecondaries(), G4DNAPTBExcitationModel::SampleSecondaries(), and G4DNAPTBIonisationModel::SampleSecondaries().

◆ ReadAndSaveCSFile()

void G4VDNAModel::ReadAndSaveCSFile ( const size_t & materialID,
const G4ParticleDefinition * p,
const G4String & file,
const G4double & scaleFactor )
protected

ReadAndSaveCSFile Read and save a "simple" cross section file : use of G4DNACrossSectionDataSet->loadData()

Parameters
materialName
particleName
file
scaleFactor

Definition at line 171 of file G4VDNAModel.cc.

173{
174 fData[materialID][p] =
175 std::make_unique<G4DNACrossSectionDataSet>(new G4LogLogInterpolation, eV, scaleFactor);
176 fData[materialID][p]->LoadData(file);
177}

Referenced by LoadCrossSectionData().

◆ ReadDiffCSFile()

void G4VDNAModel::ReadDiffCSFile ( const size_t & materialName,
const G4ParticleDefinition * particleName,
const G4String & path,
const G4double & scaleFactor )
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.

Parameters
materialName
particleName
path
scaleFactor

Definition at line 122 of file G4VDNAModel.cc.

124{
125 G4String text(
126 "ReadDiffCSFile must be implemented in the model class using a differential cross section data "
127 "file");
128
129 G4Exception("G4VDNAModel::ReadDiffCSFile", "em0003", FatalException, text);
130}

Referenced by LoadCrossSectionData().

◆ SampleSecondaries()

void G4VDNAModel::SampleSecondaries ( std::vector< G4DynamicParticle * > * ,
const G4MaterialCutsCouple * ,
const G4DynamicParticle * ,
G4double tmin = 0,
G4double tmax = DBL_MAX )
overridepure 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

param particleChangeForGamma

Parameters
tmin

param tmax

Implements G4VEmModel.

Implemented in G4DNACPA100ElasticModel, G4DNACPA100ExcitationModel, G4DNACPA100IonisationModel, G4DNAPTBElasticModel, G4DNAPTBExcitationModel, G4DNAPTBIonisationModel, and G4DNAVacuumModel.

◆ SetHighELimit()

void G4VDNAModel::SetHighELimit ( const size_t & materialID,
const G4ParticleDefinition * particle,
G4double lim )
inline

SetHighEnergyLimit.

Parameters
material
particle
lim

Definition at line 152 of file G4VDNAModel.hh.

153 {
154 fHighEnergyLimits[materialID][particle] = lim;
155 }

Referenced by G4DNACPA100ElasticModel::Initialise(), G4DNACPA100ExcitationModel::Initialise(), G4DNACPA100IonisationModel::Initialise(), G4DNAPTBElasticModel::Initialise(), G4DNAPTBExcitationModel::Initialise(), and G4DNAPTBIonisationModel::Initialise().

◆ SetLowELimit()

void G4VDNAModel::SetLowELimit ( const size_t & materialID,
const G4ParticleDefinition * particle,
G4double lim )
inline

SetLowEnergyLimit.

Parameters
material
particle
lim

Definition at line 163 of file G4VDNAModel.hh.

164 {
165 fLowEnergyLimits[materialID][particle] = lim;
166 }

Referenced by G4DNACPA100ElasticModel::Initialise(), G4DNACPA100ExcitationModel::Initialise(), G4DNACPA100IonisationModel::Initialise(), G4DNAPTBElasticModel::Initialise(), G4DNAPTBExcitationModel::Initialise(), and G4DNAPTBIonisationModel::Initialise().

Member Data Documentation

◆ isInitialised

G4bool G4VDNAModel::isInitialised = false
protected

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