Geant4 10.7.0
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=0, const G4String &nam="DNADingfelderChargeIncreaseModel")
 
virtual ~G4DNADingfelderChargeIncreaseModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double CrossSectionPerVolume (const G4Material *material, const G4ParticleDefinition *p, G4double ekin, G4double emin, G4double emax)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
void SelectStationary (G4bool input)
 
- 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 *, G4double &eloss, G4double &niel, G4double length)
 
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 ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector< G4EmElementSelector * > * GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual 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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
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 *)
 
const G4ElementGetCurrentElement () const
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

Protected Attributes

G4ParticleChangeForGammafParticleChangeForGamma
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 

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 40 of file G4DNADingfelderChargeIncreaseModel.hh.

Constructor & Destructor Documentation

◆ G4DNADingfelderChargeIncreaseModel()

G4DNADingfelderChargeIncreaseModel::G4DNADingfelderChargeIncreaseModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNADingfelderChargeIncreaseModel" 
)

Definition at line 39 of file G4DNADingfelderChargeIncreaseModel.cc.

40 :
41G4VEmModel(nam), isInitialised(false)
42{
43 fpMolWaterDensity = 0;
44
45 numberOfPartialCrossSections[0] = 0;
46 numberOfPartialCrossSections[1] = 0;
47
48 verboseLevel = 0;
49 // Verbosity scale:
50 // 0 = nothing
51 // 1 = warning for energy non-conservation
52 // 2 = details of energy budget
53 // 3 = calculation of cross sections, file openings, sampling of atoms
54 // 4 = entering in methods
55
56 if (verboseLevel > 0)
57 {
58 G4cout << "Dingfelder charge increase model is constructed " << G4endl;
59 }
61
62 // Selection of stationary mode
63
64 statCode = false;
65}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout

◆ ~G4DNADingfelderChargeIncreaseModel()

G4DNADingfelderChargeIncreaseModel::~G4DNADingfelderChargeIncreaseModel ( )
virtual

Definition at line 69 of file G4DNADingfelderChargeIncreaseModel.cc.

70{}

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 196 of file G4DNADingfelderChargeIncreaseModel.cc.

201{
202 if (verboseLevel > 3)
203 {
204 G4cout
205 << "Calling CrossSectionPerVolume() of G4DNADingfelderChargeIncreaseModel"
206 << G4endl;
207 }
208
209 // Calculate total cross section for model
210
211 G4DNAGenericIonsManager *instance;
213
214 if (
215 particleDefinition != instance->GetIon("hydrogen")
216 &&
217 particleDefinition != instance->GetIon("alpha+")
218 &&
219 particleDefinition != instance->GetIon("helium")
220 )
221
222 return 0;
223
224 G4double lowLim = 0;
225 G4double highLim = 0;
226 G4double totalCrossSection = 0.;
227
228 G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
229
230 const G4String& particleName = particleDefinition->GetParticleName();
231
232 std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
233 pos1 = lowEnergyLimit.find(particleName);
234
235 if (pos1 != lowEnergyLimit.end())
236 {
237 lowLim = pos1->second;
238 }
239
240 std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
241 pos2 = highEnergyLimit.find(particleName);
242
243 if (pos2 != highEnergyLimit.end())
244 {
245 highLim = pos2->second;
246 }
247
248 if (k >= lowLim && k <= highLim)
249 {
250 //HYDROGEN
251 if (particleDefinition == instance->GetIon("hydrogen"))
252 {
253 const G4double aa = 2.835;
254 const G4double bb = 0.310;
255 const G4double cc = 2.100;
256 const G4double dd = 0.760;
257 const G4double fac = 1.0e-18;
258 const G4double rr = 13.606 * eV;
259
260 G4double t = k / (proton_mass_c2/electron_mass_c2);
261 G4double x = t / rr;
262 G4double temp = 4.0 * pi * Bohr_radius/nm * Bohr_radius/nm * fac;
263 G4double sigmal = temp * cc * (std::pow(x,dd));
264 G4double sigmah = temp * (aa * std::log(1.0 + x) + bb) / x;
265 totalCrossSection = 1.0/(1.0/sigmal + 1.0/sigmah) *m*m;
266 }
267 else
268 {
269 totalCrossSection = Sum(k,particleDefinition);
270 }
271 }
272
273 if (verboseLevel > 2)
274 {
275 G4cout << "__________________________________" << G4endl;
276 G4cout << "G4DNADingfelderChargeIncreaseModel - XS INFO START" << G4endl;
277 G4cout << "Kinetic energy(eV)=" << k/eV << " particle : " << particleName << G4endl;
278 G4cout << "Cross section per water molecule (cm^2)=" << totalCrossSection/cm/cm << G4endl;
279 G4cout << "Cross section per water molecule (cm^-1)=" << totalCrossSection*waterDensity/(1./cm) << G4endl;
280 // G4cout << " - Cross section per water molecule (cm^-1)="
281 // << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
282 G4cout << "G4DNADingfelderChargeIncreaseModel - XS INFO END" << G4endl;
283 }
284
285 return totalCrossSection*waterDensity;
286 // return totalCrossSection*material->GetAtomicNumDensityVector()[1];
287
288}
double G4double
Definition: G4Types.hh:83
static G4DNAGenericIonsManager * Instance(void)
G4ParticleDefinition * GetIon(const G4String &name)
size_t GetIndex() const
Definition: G4Material.hh:258
const G4double pi

◆ Initialise()

void G4DNADingfelderChargeIncreaseModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 74 of file G4DNADingfelderChargeIncreaseModel.cc.

76{
77
78 if (verboseLevel > 3)
79 {
80 G4cout << "Calling G4DNADingfelderChargeIncreaseModel::Initialise()"
81 << G4endl;
82 }
83
84 // Energy limits
85
88 G4ParticleDefinition* hydrogenDef = instance->GetIon("hydrogen");
89 G4ParticleDefinition* alphaPlusDef = instance->GetIon("alpha+");
90 G4ParticleDefinition* heliumDef = instance->GetIon("helium");
91
92 G4String hydrogen;
93 G4String alphaPlus;
94 G4String helium;
95
96 // Limits
97
98 hydrogen = hydrogenDef->GetParticleName();
99 lowEnergyLimit[hydrogen] = 100. * eV;
100 highEnergyLimit[hydrogen] = 100. * MeV;
101
102 alphaPlus = alphaPlusDef->GetParticleName();
103 lowEnergyLimit[alphaPlus] = 1. * keV;
104 highEnergyLimit[alphaPlus] = 400. * MeV;
105
106 helium = heliumDef->GetParticleName();
107 lowEnergyLimit[helium] = 1. * keV;
108 highEnergyLimit[helium] = 400. * MeV;
109
110 //
111
112 if (particle==hydrogenDef)
113 {
114 SetLowEnergyLimit(lowEnergyLimit[hydrogen]);
115 SetHighEnergyLimit(highEnergyLimit[hydrogen]);
116 }
117
118 if (particle==alphaPlusDef)
119 {
120 SetLowEnergyLimit(lowEnergyLimit[alphaPlus]);
121 SetHighEnergyLimit(highEnergyLimit[alphaPlus]);
122 }
123
124 if (particle==heliumDef)
125 {
126 SetLowEnergyLimit(lowEnergyLimit[helium]);
127 SetHighEnergyLimit(highEnergyLimit[helium]);
128 }
129
130 // Final state
131
132 //ALPHA+
133
134 f0[0][0]=1.;
135 a0[0][0]=2.25;
136 a1[0][0]=-0.75;
137 b0[0][0]=-32.10;
138 c0[0][0]=0.600;
139 d0[0][0]=2.40;
140 x0[0][0]=4.60;
141
142 x1[0][0]=-1.;
143 b1[0][0]=-1.;
144
145 numberOfPartialCrossSections[0]=1;
146
147 //HELIUM
148
149 f0[0][1]=1.;
150 a0[0][1]=2.25;
151 a1[0][1]=-0.75;
152 b0[0][1]=-30.93;
153 c0[0][1]=0.590;
154 d0[0][1]=2.35;
155 x0[0][1]=4.29;
156
157 f0[1][1]=1.;
158 a0[1][1]=2.25;
159 a1[1][1]=-0.75;
160 b0[1][1]=-32.61;
161 c0[1][1]=0.435;
162 d0[1][1]=2.70;
163 x0[1][1]=4.45;
164
165 x1[0][1]=-1.;
166 b1[0][1]=-1.;
167
168 x1[1][1]=-1.;
169 b1[1][1]=-1.;
170
171 numberOfPartialCrossSections[1]=2;
172
173 //
174
175 if( verboseLevel>0 )
176 {
177 G4cout << "Dingfelder charge increase model is initialized " << G4endl
178 << "Energy range: "
179 << LowEnergyLimit() / keV << " keV - "
180 << HighEnergyLimit() / MeV << " MeV for "
181 << particle->GetParticleName()
182 << G4endl;
183 }
184
185 // Initialize water density pointer
187
188 if (isInitialised) return;
189
191 isInitialised = true;
192}
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)
Definition: G4Material.cc:651
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:764

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 292 of file G4DNADingfelderChargeIncreaseModel.cc.

298{
299 if (verboseLevel > 3)
300 {
301 G4cout
302 << "Calling SampleSecondaries() of G4DNADingfelderChargeIncreaseModel"
303 << G4endl;
304 }
305
307
308 G4ParticleDefinition* definition = aDynamicParticle->GetDefinition();
309
310 G4double particleMass = definition->GetPDGMass();
311
312 G4double inK = aDynamicParticle->GetKineticEnergy();
313
314 G4int finalStateIndex = RandomSelect(inK,definition);
315
316 G4int n = NumberOfFinalStates(definition,finalStateIndex);
317
318 G4double outK = 0.;
319
320 if (!statCode) outK = inK - IncomingParticleBindingEnergyConstant(definition,finalStateIndex);
321
322 else outK = inK;
323
324 if (statCode)
326 ProposeLocalEnergyDeposit(IncomingParticleBindingEnergyConstant(definition,finalStateIndex));
327
329
330 G4DNAGenericIonsManager* instance;
332
333 G4double electronK;
334 if (definition == instance->GetIon("hydrogen")) electronK = inK*electron_mass_c2/proton_mass_c2;
335 else electronK = inK*electron_mass_c2/(particleMass);
336
337 if (outK<0)
338 {
339 G4Exception("G4DNADingfelderChargeIncreaseModel::SampleSecondaries","em0004",
340 FatalException,"Final kinetic energy is negative.");
341 }
342
343 G4DynamicParticle* dp = new G4DynamicParticle(OutgoingParticleDefinition(definition,finalStateIndex),
344 aDynamicParticle->GetMomentumDirection(),
345 outK);
346
347 fvect->push_back(dp);
348
349 n = n - 1;
350
351 while (n>0)
352 {
353 n--;
354 fvect->push_back(new G4DynamicParticle
355 (G4Electron::Electron(), aDynamicParticle->GetMomentumDirection(), electronK) );
356 }
357}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
@ fStopAndKill
int G4int
Definition: G4Types.hh:85
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition: G4Electron.cc:93
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

◆ SelectStationary()

void G4DNADingfelderChargeIncreaseModel::SelectStationary ( G4bool  input)
inline

Definition at line 121 of file G4DNADingfelderChargeIncreaseModel.hh.

122{
123 statCode = input;
124}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNADingfelderChargeIncreaseModel::fParticleChangeForGamma
protected

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