Geant4 11.2.2
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.
 
 ~G4DNAModelInterface () override=default
 ~G4DNAModelManager Destructor
 
 G4DNAModelInterface (const G4DNAModelInterface &)=delete
 
G4DNAModelInterfaceoperator= (const G4DNAModelInterface &right)=delete
 
void Initialise (const G4ParticleDefinition *particle, const G4DataVector &cuts) override
 Initialise Initialise method to call all the initialise methods of the registered models.
 
G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
 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).
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *fVect, const G4MaterialCutsCouple *couple, const G4DynamicParticle *aDynamicElectron, G4double tmin, G4double tmax) override
 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 (G4VEmModel *model)
 RegisterModel Method used to associate a model with the interaction.
 
std::size_t GetSelectedMaterial ()
 GetSelectedMaterial To allow the user to retrieve the selected material in case of a composite material.
 
void StreamInfo (std::ostream &os) const
 
- 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
 

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
 
std::size_t currentCoupleIndex = 0
 
std::size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

Detailed Description

Definition at line 40 of file G4DNAModelInterface.hh.

Constructor & Destructor Documentation

◆ G4DNAModelInterface() [1/2]

G4DNAModelInterface::G4DNAModelInterface ( const G4String & nam)
explicit

G4DNAModelManager Constructor.

Parameters
nam

Definition at line 38 of file G4DNAModelInterface.cc.

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

◆ ~G4DNAModelInterface()

G4DNAModelInterface::~G4DNAModelInterface ( )
overridedefault

~G4DNAModelManager Destructor

◆ G4DNAModelInterface() [2/2]

G4DNAModelInterface::G4DNAModelInterface ( const G4DNAModelInterface & )
delete

Member Function Documentation

◆ CrossSectionPerVolume()

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

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 73 of file G4DNAModelInterface.cc.

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

◆ GetSelectedMaterial()

std::size_t 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 105 of file G4DNAModelInterface.hh.

106 {
107 return fSampledMat;
108 }

◆ Initialise()

void G4DNAModelInterface::Initialise ( const G4ParticleDefinition * particle,
const G4DataVector & cuts )
overridevirtual

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

Parameters
particle
cuts

Implements G4VEmModel.

Definition at line 42 of file G4DNAModelInterface.cc.

43{
44 // Those two statements are necessary to override the energy limits set in the G4DNAProcesses
45 // (ionisation, elastic, etc...). Indeed, with the ModelInterface system, the model define
46 // themselves their energy limits per material and particle. Therefore, such a limit should not be
47 // in the G4DNAProcess classes.
48 //
49
50 fpG4_WATER = G4Material::GetMaterial("G4_WATER", false);
51
54
55 fpParticleChangeForGamma = GetParticleChangeForGamma();
56
57 // Loop on all the registered models to initialise them
58 for (auto & fRegisteredModel : fRegisteredModels) {
59 fRegisteredModel->SetParticleChange(fpParticleChangeForGamma);
60 fRegisteredModel->Initialise(particle, cuts);
61 }
62 // used to retrieve the model corresponding to the current material/particle couple
63 BuildMaterialParticleModelTable(particle);
64
65 BuildMaterialMolPerVolTable();
66
68
69}
void StreamInfo(std::ostream &os) const
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
void SetHighEnergyLimit(G4double)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
void SetLowEnergyLimit(G4double)
#define DBL_MAX
Definition templates.hh:62

◆ operator=()

G4DNAModelInterface & G4DNAModelInterface::operator= ( const G4DNAModelInterface & right)
delete

◆ RegisterModel()

void G4DNAModelInterface::RegisterModel ( G4VEmModel * model)

RegisterModel Method used to associate a model with the interaction.

Parameters
model

Definition at line 306 of file G4DNAModelInterface.cc.

307{
308 fRegisteredModels.push_back(model);
309}

◆ SampleSecondaries()

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

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

param couple

Parameters
aDynamicElectron

param tmin

Parameters
tmax

Implements G4VEmModel.

Definition at line 206 of file G4DNAModelInterface.cc.

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

◆ StreamInfo()

void G4DNAModelInterface::StreamInfo ( std::ostream & os) const

Definition at line 499 of file G4DNAModelInterface.cc.

500{
501 G4long prec = os.precision(5);
502 os << "======================================= Materials of " << std::setw(17)
503 << this->GetName() << " ================================================"
504 << "\n";
505 os << std::setw(15) << "Material#" << std::setw(13) << "Particle" << std::setw(35) << "Model"
506 << std::setw(17) << "LowLimit(MeV)" << std::setw(17) << "HighLimit(MeV)" << std::setw(13)
507 << "Fast" << std::setw(13) << "Stationary" << std::setw(13) << "Chemistry" << G4endl;
508 for (const auto& it1 : fMaterialParticleModelTable) {
509 os << std::setw(15) << (*G4Material::GetMaterialTable())[it1.first]->GetName();
510 for (const auto& it2 : it1.second) {
511 os << std::setw(13) << it2.first->GetParticleName();
512 os << std::setw(35) << it2.second->GetName();
513 auto DNAModel = dynamic_cast<G4VDNAModel*>(it2.second);
514 if (DNAModel == nullptr) {
515 os << std::setw(17) << it2.second->LowEnergyLimit();
516 os << std::setw(17) << it2.second->HighEnergyLimit();
517 }
518 else {
519 auto lowL = DNAModel->GetLowELimit(it1.first, it2.first);
520 auto highL = DNAModel->GetHighELimit(it1.first, it2.first);
521 os << std::setw(17) << lowL;
522 os << std::setw(17) << highL;
523 }
524 os << std::setw(13) << "no";
525 os << std::setw(13) << "no";
526 os << std::setw(13) << "no" << G4endl;
527 }
528 }
529
530 os << "========================================================================================"
531 "=================================================="
532 << G4endl;
533 os.precision(prec);
534}
long G4long
Definition G4Types.hh:87
static G4MaterialTable * GetMaterialTable()
const G4String & GetName() const

Referenced by Initialise().


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