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

#include <G4DNADingfelderChargeIncreaseModel.hh>

+ Inheritance diagram for G4DNADingfelderChargeIncreaseModel:

Public Member Functions

 G4DNADingfelderChargeIncreaseModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="DNADingfelderChargeIncreaseModel")
 
 ~G4DNADingfelderChargeIncreaseModel () override=default
 
G4DNADingfelderChargeIncreaseModeloperator= (const G4DNADingfelderChargeIncreaseModel &right)=delete
 
 G4DNADingfelderChargeIncreaseModel (const G4DNADingfelderChargeIncreaseModel &)=delete
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
void SelectStationary (G4bool input)
 
- 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 Attributes

G4ParticleChangeForGammafParticleChangeForGamma = nullptr
 
- 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
 

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 *)
 

Detailed Description

Definition at line 39 of file G4DNADingfelderChargeIncreaseModel.hh.

Constructor & Destructor Documentation

◆ G4DNADingfelderChargeIncreaseModel() [1/2]

G4DNADingfelderChargeIncreaseModel::G4DNADingfelderChargeIncreaseModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "DNADingfelderChargeIncreaseModel" )
explicit

Definition at line 44 of file G4DNADingfelderChargeIncreaseModel.cc.

45 :
46G4VEmModel(nam)
47{
48 if (verboseLevel > 0)
49 {
50 G4cout << "Dingfelder charge increase model is constructed " << G4endl;
51 }
52}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67

◆ ~G4DNADingfelderChargeIncreaseModel()

G4DNADingfelderChargeIncreaseModel::~G4DNADingfelderChargeIncreaseModel ( )
overridedefault

◆ G4DNADingfelderChargeIncreaseModel() [2/2]

G4DNADingfelderChargeIncreaseModel::G4DNADingfelderChargeIncreaseModel ( const G4DNADingfelderChargeIncreaseModel & )
delete

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 179 of file G4DNADingfelderChargeIncreaseModel.cc.

184{
185 if (verboseLevel > 3)
186 {
187 G4cout
188 << "Calling CrossSectionPerVolume() of G4DNADingfelderChargeIncreaseModel"
189 << G4endl;
190 }
191
192 // Calculate total cross section for model
193
194 if (
195 particleDefinition != hydrogenDef
196 &&
197 particleDefinition != alphaPlusDef
198 &&
199 particleDefinition != heliumDef
200 )
201
202 return 0;
203
204 G4double lowLim = 0;
205 G4double highLim = 0;
206 G4double totalCrossSection = 0.;
207
208 G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
209
210 const G4String& particleName = particleDefinition->GetParticleName();
211
212 std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
213 pos1 = lowEnergyLimit.find(particleName);
214
215 if (pos1 != lowEnergyLimit.end())
216 {
217 lowLim = pos1->second;
218 }
219
220 std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
221 pos2 = highEnergyLimit.find(particleName);
222
223 if (pos2 != highEnergyLimit.end())
224 {
225 highLim = pos2->second;
226 }
227
228 if (k >= lowLim && k <= highLim)
229 {
230 //HYDROGEN
231 if (particleDefinition == hydrogenDef)
232 {
233 const G4double aa = 2.835;
234 const G4double bb = 0.310;
235 const G4double cc = 2.100;
236 const G4double dd = 0.760;
237 const G4double fac = 1.0e-18;
238 const G4double rr = 13.606 * eV;
239
240 G4double t = k / (proton_mass_c2/electron_mass_c2);
241 G4double x = t / rr;
242 G4double temp = 4.0 * pi * Bohr_radius/nm * Bohr_radius/nm * fac;
243 G4double sigmal = temp * cc * (gpow->powA(x,dd));
244 G4double sigmah = temp * (aa * G4Log(1.0 + x) + bb) / x;
245 totalCrossSection = 1.0/(1.0/sigmal + 1.0/sigmah) *m*m;
246 }
247 else
248 {
249 totalCrossSection = Sum(k,particleDefinition);
250 }
251 }
252
253 if (verboseLevel > 2)
254 {
255 G4cout << "__________________________________" << G4endl;
256 G4cout << "G4DNADingfelderChargeIncreaseModel - XS INFO START" << G4endl;
257 G4cout << "Kinetic energy(eV)=" << k/eV << " particle : " << particleName << G4endl;
258 G4cout << "Cross section per water molecule (cm^2)=" << totalCrossSection/cm/cm << G4endl;
259 G4cout << "Cross section per water molecule (cm^-1)=" << totalCrossSection*waterDensity/(1./cm) << G4endl;
260 // G4cout << " - Cross section per water molecule (cm^-1)="
261 // << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
262 G4cout << "G4DNADingfelderChargeIncreaseModel - XS INFO END" << G4endl;
263 }
264
265 return totalCrossSection*waterDensity;
266 // return totalCrossSection*material->GetAtomicNumDensityVector()[1];
267
268}
G4double G4Log(G4double x)
Definition G4Log.hh:227
double G4double
Definition G4Types.hh:83
std::size_t GetIndex() const
G4double powA(G4double A, G4double y) const
Definition G4Pow.hh:230

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 56 of file G4DNADingfelderChargeIncreaseModel.cc.

58{
59
60 if (verboseLevel > 3)
61 {
62 G4cout << "Calling G4DNADingfelderChargeIncreaseModel::Initialise()"
63 << G4endl;
64 }
65
66 // Energy limits
67
70 hydrogenDef = instance->GetIon("hydrogen");
71 alphaPlusPlusDef = G4Alpha::Alpha();
72 alphaPlusDef = instance->GetIon("alpha+");
73 heliumDef = instance->GetIon("helium");
74
75 G4String hydrogen;
76 G4String alphaPlus;
77 G4String helium;
78
79 // Limits
80
81 hydrogen = hydrogenDef->GetParticleName();
82 lowEnergyLimit[hydrogen] = 100. * eV;
83 highEnergyLimit[hydrogen] = 100. * MeV;
84
85 alphaPlus = alphaPlusDef->GetParticleName();
86 lowEnergyLimit[alphaPlus] = 1. * keV;
87 highEnergyLimit[alphaPlus] = 400. * MeV;
88
89 helium = heliumDef->GetParticleName();
90 lowEnergyLimit[helium] = 1. * keV;
91 highEnergyLimit[helium] = 400. * MeV;
92
93 //
94
95 if (particle==hydrogenDef)
96 {
97 SetLowEnergyLimit(lowEnergyLimit[hydrogen]);
98 SetHighEnergyLimit(highEnergyLimit[hydrogen]);
99 }
100
101 if (particle==alphaPlusDef)
102 {
103 SetLowEnergyLimit(lowEnergyLimit[alphaPlus]);
104 SetHighEnergyLimit(highEnergyLimit[alphaPlus]);
105 }
106
107 if (particle==heliumDef)
108 {
109 SetLowEnergyLimit(lowEnergyLimit[helium]);
110 SetHighEnergyLimit(highEnergyLimit[helium]);
111 }
112
113 // Final state
114
115 //ALPHA+
116
117 f0[0][0]=1.;
118 a0[0][0]=2.25;
119 a1[0][0]=-0.75;
120 b0[0][0]=-32.10;
121 c0[0][0]=0.600;
122 d0[0][0]=2.40;
123 x0[0][0]=4.60;
124
125 x1[0][0]=-1.;
126 b1[0][0]=-1.;
127
128 numberOfPartialCrossSections[0]=1;
129
130 //HELIUM
131
132 f0[0][1]=1.;
133 a0[0][1]=2.25;
134 a1[0][1]=-0.75;
135 b0[0][1]=-30.93;
136 c0[0][1]=0.590;
137 d0[0][1]=2.35;
138 x0[0][1]=4.29;
139
140 f0[1][1]=1.;
141 a0[1][1]=2.25;
142 a1[1][1]=-0.75;
143 b0[1][1]=-32.61;
144 c0[1][1]=0.435;
145 d0[1][1]=2.70;
146 x0[1][1]=4.45;
147
148 x1[0][1]=-1.;
149 b1[0][1]=-1.;
150
151 x1[1][1]=-1.;
152 b1[1][1]=-1.;
153
154 numberOfPartialCrossSections[1]=2;
155
156 //
157
158 if( verboseLevel>0 )
159 {
160 G4cout << "Dingfelder charge increase model is initialized " << G4endl
161 << "Energy range: "
162 << LowEnergyLimit() / keV << " keV - "
163 << HighEnergyLimit() / MeV << " MeV for "
164 << particle->GetParticleName()
165 << G4endl;
166 }
167
168 // Initialize water density pointer
170
171 if (isInitialised) return;
172
174 isInitialised = true;
175}
static G4Alpha * Alpha()
Definition G4Alpha.cc:83
static G4DNAGenericIonsManager * Instance()
G4ParticleDefinition * GetIon(const G4String &name)
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()
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
G4double LowEnergyLimit() const
G4double HighEnergyLimit() const
void SetLowEnergyLimit(G4double)

◆ operator=()

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

◆ SampleSecondaries()

void G4DNADingfelderChargeIncreaseModel::SampleSecondaries ( std::vector< G4DynamicParticle * > * fvect,
const G4MaterialCutsCouple * ,
const G4DynamicParticle * aDynamicParticle,
G4double tmin,
G4double maxEnergy )
overridevirtual

Implements G4VEmModel.

Definition at line 272 of file G4DNADingfelderChargeIncreaseModel.cc.

278{
279 if (verboseLevel > 3)
280 {
281 G4cout
282 << "Calling SampleSecondaries() of G4DNADingfelderChargeIncreaseModel"
283 << G4endl;
284 }
285
287
288 G4ParticleDefinition* definition = aDynamicParticle->GetDefinition();
289
290 G4double particleMass = definition->GetPDGMass();
291
292 G4double inK = aDynamicParticle->GetKineticEnergy();
293
294 G4int finalStateIndex = RandomSelect(inK,definition);
295
296 G4int n = NumberOfFinalStates(definition,finalStateIndex);
297
298 G4double outK = 0.;
299
300 if (!statCode) outK = inK - IncomingParticleBindingEnergyConstant(definition,finalStateIndex);
301
302 else outK = inK;
303
304 if (statCode)
306 ProposeLocalEnergyDeposit(IncomingParticleBindingEnergyConstant(definition,finalStateIndex));
307
309
310 G4double electronK;
311 if (definition == hydrogenDef) electronK = inK*electron_mass_c2/proton_mass_c2;
312 else electronK = inK*electron_mass_c2/(particleMass);
313
314 if (outK<0)
315 {
316 G4Exception("G4DNADingfelderChargeIncreaseModel::SampleSecondaries","em0004",
317 FatalException,"Final kinetic energy is negative.");
318 }
319
320 auto dp = new G4DynamicParticle(OutgoingParticleDefinition(definition,finalStateIndex),
321 aDynamicParticle->GetMomentumDirection(),
322 outK);
323
324 fvect->push_back(dp);
325
326 n = n - 1;
327
328 while (n>0)
329 {
330 n--;
331 fvect->push_back(new G4DynamicParticle
332 (G4Electron::Electron(), aDynamicParticle->GetMomentumDirection(), electronK) );
333 }
334}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
@ fStopAndKill
int G4int
Definition G4Types.hh:85
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition G4Electron.cc:91
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

◆ SelectStationary()

void G4DNADingfelderChargeIncreaseModel::SelectStationary ( G4bool input)
inline

Definition at line 112 of file G4DNADingfelderChargeIncreaseModel.hh.

113{
114 statCode = input;
115}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNADingfelderChargeIncreaseModel::fParticleChangeForGamma = nullptr
protected

Definition at line 63 of file G4DNADingfelderChargeIncreaseModel.hh.

Referenced by Initialise(), and SampleSecondaries().


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