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

#include <G4DNAModelInterface.hh>

+ Inheritance diagram for G4DNAModelInterface:

Public Member Functions

 G4DNAModelInterface (const G4String &nam)
 G4DNAModelManager Constructor.
 
virtual ~G4DNAModelInterface ()
 ~G4DNAModelManager Destructor
 
virtual void Initialise (const G4ParticleDefinition *particle, const G4DataVector &cuts)
 Initialise Initialise method to call all the initialise methods of the registered models.
 
virtual G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
 CrossSectionPerVolume Method called by the process and used to call the CrossSectionPerVolume method of the registered models. The method also calculates through G4DNAMolecularMaterial the number of molecule per volume unit for the current material or (component of a composite material).
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *fVect, const G4MaterialCutsCouple *couple, const G4DynamicParticle *aDynamicElectron, G4double tmin, G4double tmax)
 SampleSecondaries Used to call the SampleSecondaries method of the registered models. A sampling is done to select a component if the material is a composite one.
 
void RegisterModel (G4VDNAModel *model)
 RegisterModel Method used to associate a model with the interaction.
 
void RegisterModel (G4VEmModel *model, const G4ParticleDefinition *particle)
 
G4String GetSelectedMaterial ()
 GetSelectedMaterial To allow the user to retrieve the selected material in case of a composite material.
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
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 CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=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 LPMFlag () 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 SetLPMFlag (G4bool val)
 
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)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 
- 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
 
size_t currentCoupleIndex = 0
 
size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

Detailed Description

Definition at line 45 of file G4DNAModelInterface.hh.

Constructor & Destructor Documentation

◆ G4DNAModelInterface()

G4DNAModelInterface::G4DNAModelInterface ( const G4String nam)

G4DNAModelManager Constructor.

Parameters
nam

Definition at line 37 of file G4DNAModelInterface.cc.

38 : G4VEmModel(nam), fName(nam), fpParticleChangeForGamma(0), fSampledMat("")
39{
40
41}

◆ ~G4DNAModelInterface()

G4DNAModelInterface::~G4DNAModelInterface ( )
virtual

~G4DNAModelManager Destructor

Definition at line 45 of file G4DNAModelInterface.cc.

46{
47 // Loop on all the registered models to properly delete them (free the memory)
48 for(std::size_t i=0, ie = fRegisteredModels.size(); i<ie; ++i)
49 {
50 if(fRegisteredModels.at(i) != nullptr) delete fRegisteredModels.at(i);
51 }
52}

Member Function Documentation

◆ CrossSectionPerVolume()

G4double G4DNAModelInterface::CrossSectionPerVolume ( const G4Material material,
const G4ParticleDefinition p,
G4double  ekin,
G4double  emin,
G4double  emax 
)
virtual

CrossSectionPerVolume Method called by the process and used to call the CrossSectionPerVolume method of the registered models. The method also calculates through G4DNAMolecularMaterial the number of molecule per volume unit for the current material or (component of a composite material).

Parameters
material
p
ekin
emin
emax
Returns
the final cross section value times with the number of molecule per volume unit

Reimplemented from G4VEmModel.

Definition at line 84 of file G4DNAModelInterface.cc.

89{
90 // Method to return the crossSection * nbMoleculePerUnitVolume to the process class.
91 // Process class then calculates the path.
92 // The cross section is calculated in the registered model(s) and this class just call the method
93 // Two cases are handled here: normal material and composite material.
94 //
95 // Idea:
96 // *** Simple material ***
97 // Ask for the cross section of the chosen model.
98 // Multiply it by the number of medium molecules per volume unit.
99 // Return the value.
100 // *** Composite material ***
101 // Ask for the cross section of the chosen model for each component.
102 // Apply a factor to each cross section and sum the results. The factor is the molecule number of component per composite volume unit.
103 // The total cross section is returned.
104
105 // To reset the sampledMat variable.
106 // Can be used by user to retrieve current component
107 fSampledMat = "";
108
109 // This is the value to be sum up and to be returned at then end
110 G4double crossSectionTimesNbMolPerVol (0);
111
112 // Reset the map saving the material and the cumulated corresponding cross section
113 // Used in SampleSecondaries if the interaction is selected for the step and if the material is a composite
114 fMaterialCS.clear();
115
116 // This is the value to be used by SampleSecondaries
117 fCSsumTot = 0;
118
119 // *****************************
120 // Material is not a composite
121 // *****************************
122 //
123 if(material->GetMatComponents().empty())
124 {
125 // Get the material name
126 const G4String& materialName = material->GetName();
127
128 // Use the table to get the model
129 G4VDNAModel* model = GetDNAModel(materialName, p->GetParticleName(), ekin);
130
131 // Get the nunber of molecules per volume unit for that material
132 G4double nbOfMoleculePerVolumeUnit = GetNumMoleculePerVolumeUnitForMaterial(material);
133
134 // Calculate the cross section times the number of molecules
135 if(model != 0)
136 crossSectionTimesNbMolPerVol = nbOfMoleculePerVolumeUnit * model->CrossSectionPerVolume(material, materialName, p, ekin, emin, emax);
137 else // no model was selected, we are out of the energy ranges
138 crossSectionTimesNbMolPerVol = 0.;
139 }
140
141 // ********************************
142 // Material is a composite
143 // ********************************
144 //
145 else
146 {
147 // Copy the map in a local variable
148 // Otherwise we get segmentation fault and iterator pointing to nowhere: do not know why...
149 // Maybe MatComponents map is overrided by something somewhere ?
150 std::map<G4Material*, G4double> componentsMap = material->GetMatComponents();
151
152 // Retrieve the iterator
153 std::map<G4Material*, G4double>::const_iterator it = componentsMap.begin();
154
155 // Get the size
156 std::size_t componentNumber = componentsMap.size();
157
158 // Loop on all the components
159 //for(it = material->GetMatComponents().begin(); it!=material->GetMatComponents().end();++it)
160 for(std::size_t i=0; i<componentNumber; ++i)
161 {
162 // Get the current component
163 G4Material* component = it->first;
164
165 // Get the current component mass fraction
166 //G4double massFraction = it->second;
167
168 // Get the number of component molecules in a volume unit of composite material
169 G4double nbMoleculeOfComponentInCompositeMat = GetNumMolPerVolUnitForComponentInComposite(component, material);
170
171 // Get the current component name
172 const G4String componentName = component->GetName();
173
174 // Retrieve the model corresponding to the current component (ie material)
175 G4VDNAModel* model = GetDNAModel(componentName, p->GetParticleName(), ekin);
176
177 // Add the component part of the cross section to the cross section variable.
178 // The component cross section is multiplied by the total molecule number in the composite scaled by the mass fraction.
179 if(model != 0)
180 crossSectionTimesNbMolPerVol =
181 nbMoleculeOfComponentInCompositeMat * model->CrossSectionPerVolume(component, componentName, p, ekin, emin, emax);
182 else // no model was selected, we are out of the energy ranges
183 crossSectionTimesNbMolPerVol = 0.;
184
185 // Save the component name and its calculated crossSectionTimesNbMolPerVol
186 // To be used by sampling secondaries if the interaction is selected for the step
187 fMaterialCS[componentName] = crossSectionTimesNbMolPerVol;
188
189 // Save the component name and its calculated crossSectionTimesNbMolPerVol
190 // To be used by sampling secondaries if the interaction is selected for the step
191 fCSsumTot += crossSectionTimesNbMolPerVol;
192
193 // Move forward the iterator
194 ++it;
195 }
196
197 crossSectionTimesNbMolPerVol = fCSsumTot;
198
199 }
200
201 // return the cross section times the number of molecules
202 // the path of the interaction will be calculated using that value
203 return crossSectionTimesNbMolPerVol;
204}
double G4double
Definition: G4Types.hh:83
const std::map< G4Material *, G4double > & GetMatComponents() const
Definition: G4Material.hh:232
const G4String & GetName() const
Definition: G4Material.hh:172
const G4String & GetParticleName() const
The G4VDNAModel class.
Definition: G4VDNAModel.hh:50
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....

◆ GetSelectedMaterial()

G4String G4DNAModelInterface::GetSelectedMaterial ( )
inline

GetSelectedMaterial To allow the user to retrieve the selected material in case of a composite material.

Returns
the last selected material by SampleSecondaries.

Definition at line 119 of file G4DNAModelInterface.hh.

119{return fSampledMat;}

◆ Initialise()

void G4DNAModelInterface::Initialise ( const G4ParticleDefinition particle,
const G4DataVector cuts 
)
virtual

Initialise Initialise method to call all the initialise methods of the registered models.

Parameters
particle
cuts

Implements G4VEmModel.

Definition at line 56 of file G4DNAModelInterface.cc.

58{
59 // Those two statements are necessary to override the energy limits set in the G4DNAProcesses (ionisation, elastic, etc...).
60 // Indeed, with the ModelInterface system, the model define themselves their energy limits per material and particle.
61 // Therefore, such a limit should not be in the G4DNAProcess classes.
62 //
65
66 fpParticleChangeForGamma = GetParticleChangeForGamma();
67
68 // Loop on all the registered models to initialise them
69 for(std::size_t i=0, ie = fRegisteredModels.size(); i<ie; ++i)
70 {
71 fRegisteredModels.at(i)->Initialise(particle, cuts, fpParticleChangeForGamma);
72 }
73
74
75 // Build the [material][particle]=Models table
76 // used to retrieve the model corresponding to the current material/particle couple
77 BuildMaterialParticleModelTable(particle);
78
79 BuildMaterialMolPerVolTable();
80}
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:124
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:753
#define DBL_MAX
Definition: templates.hh:62

◆ RegisterModel() [1/2]

void G4DNAModelInterface::RegisterModel ( G4VDNAModel model)

RegisterModel Method used to associate a model with the interaction.

Parameters
model

Definition at line 316 of file G4DNAModelInterface.cc.

317{
318 fRegisteredModels.push_back(model);
319}

Referenced by RegisterModel().

◆ RegisterModel() [2/2]

void G4DNAModelInterface::RegisterModel ( G4VEmModel model,
const G4ParticleDefinition particle 
)

Definition at line 321 of file G4DNAModelInterface.cc.

322{
323 G4DNADummyModel* dummyWrapper = new G4DNADummyModel("G4_WATER", particle, model->GetName(), model);
324
325 RegisterModel(dummyWrapper);
326}
void RegisterModel(G4VDNAModel *model)
RegisterModel Method used to associate a model with the interaction.
const G4String & GetName() const
Definition: G4VEmModel.hh:816

◆ SampleSecondaries()

void G4DNAModelInterface::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fVect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle aDynamicElectron,
G4double  tmin,
G4double  tmax 
)
virtual

SampleSecondaries Used to call the SampleSecondaries method of the registered models. A sampling is done to select a component if the material is a composite one.

Parameters
fVect
couple
aDynamicElectron
tmin
tmax

Implements G4VEmModel.

Definition at line 208 of file G4DNAModelInterface.cc.

213{
214 // To call the sampleSecondaries method of the registered model(s)
215 // In the case of composite material, we need to choose a component to call the method from.
216 // To do so we use a random sampling on the crossSectionTimesNbMolPerVol used in CrossSectionPerVolume method.
217 // If we enter that method it means the corresponding interaction (and process) has been chosen for the current step.
218
219 G4String materialName;
220
221 // *******************************
222 // Material is not a composite
223 // *******************************
224 //
225 if(couple->GetMaterial()->GetMatComponents().empty())
226 {
227 materialName = couple->GetMaterial()->GetName();
228 }
229
230 // ****************************
231 // Material is a composite
232 // ****************************
233 //
234 else
235 {
236 // Material is a composite
237 // We need to select a component
238
239 // We select a random number between 0 and fCSSumTot
240 G4double rand = G4UniformRand()*fCSsumTot;
241
242 G4double cumulCS (0);
243
244 G4bool result = false;
245
246 // We loop on each component cumulated cross section
247 //
248 // Retrieve the iterators
249 std::map<const G4String , G4double>::const_iterator it = fMaterialCS.begin();
250 std::map<const G4String , G4double>::const_iterator ite = fMaterialCS.end();
251 // While this is true we do not have found our component.
252 while(rand>cumulCS)
253 {
254 // Check if the sampling is ok
255 if(it==ite)
256 {
257 G4Exception("G4DNAModelManager::SampleSecondaries","em0006",
259 "The random component selection has failed: we ran into the end of the map without having a selected component");
260 return; // to make some compilers happy
261 }
262
263 // Set the cumulated value for the iteration
264 cumulCS += it->second;
265
266 // Check if we have reach the material to be selected
267 // The DBL_MAX is here to take into account a return DBL_MAX in CSPerVol for the elastic model
268 // to force elastic sampleSecondaries where the particle can be killed.
269 // Used when paticle energy is lower than limit.
270 if(rand<cumulCS || cumulCS >= DBL_MAX)
271 {
272 // we have our selected material
273 materialName = it->first;
274 result = true;
275 break;
276 }
277
278 // make the iterator move forward
279 ++it;
280 }
281
282 // Check that we get a result
283 if(!result)
284 {
285 // it is possible to end up here if the return DBL_MAX of CSPerVol in the elastic model is not taken into account
286
287 G4Exception("G4DNAModelManager::SampleSecondaries","em0006",
289 "The random component selection has failed: while loop ended without a selected component.");
290 return; // to make some compilers happy
291 }
292
293 }
294
295 // **************************************
296 // Call the SampleSecondaries method
297 // **************************************
298
299 // Rename material if modified NIST material
300 // This is needed when material is obtained from G4MaterialCutsCouple
301 if(materialName.find("_MODIFIED")!=G4String::npos)
302 {
303 materialName = materialName.substr(0,materialName.size()-9);
304 }
305
306 fSampledMat = materialName;
307
308 G4VDNAModel* model = GetDNAModel(materialName,
309 aDynamicParticle->GetParticleDefinition()->GetParticleName(),
310 aDynamicParticle->GetKineticEnergy() );
311 //fMaterialParticleModelTable[materialName][aDynamicParticle->GetDefinition()->GetParticleName()][0];
312
313 model->SampleSecondaries(fVect, couple, materialName, aDynamicParticle, fpParticleChangeForGamma, tmin, tmax);
314}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:59
bool G4bool
Definition: G4Types.hh:86
#define G4UniformRand()
Definition: Randomize.hh:52
const G4Material * GetMaterial() const
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...

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