Geant4 11.1.1
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 *, 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
 

Protected Attributes

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

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 44 of file G4DNADingfelderChargeIncreaseModel.cc.

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

◆ ~G4DNADingfelderChargeIncreaseModel()

G4DNADingfelderChargeIncreaseModel::~G4DNADingfelderChargeIncreaseModel ( )
virtual

Definition at line 74 of file G4DNADingfelderChargeIncreaseModel.cc.

75{}

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 202 of file G4DNADingfelderChargeIncreaseModel.cc.

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

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 79 of file G4DNADingfelderChargeIncreaseModel.cc.

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 295 of file G4DNADingfelderChargeIncreaseModel.cc.

301{
302 if (verboseLevel > 3)
303 {
304 G4cout
305 << "Calling SampleSecondaries() of G4DNADingfelderChargeIncreaseModel"
306 << G4endl;
307 }
308
310
311 G4ParticleDefinition* definition = aDynamicParticle->GetDefinition();
312
313 G4double particleMass = definition->GetPDGMass();
314
315 G4double inK = aDynamicParticle->GetKineticEnergy();
316
317 G4int finalStateIndex = RandomSelect(inK,definition);
318
319 G4int n = NumberOfFinalStates(definition,finalStateIndex);
320
321 G4double outK = 0.;
322
323 if (!statCode) outK = inK - IncomingParticleBindingEnergyConstant(definition,finalStateIndex);
324
325 else outK = inK;
326
327 if (statCode)
329 ProposeLocalEnergyDeposit(IncomingParticleBindingEnergyConstant(definition,finalStateIndex));
330
332
333 G4double electronK;
334 if (definition == hydrogenDef) 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:59
@ 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 127 of file G4DNADingfelderChargeIncreaseModel.hh.

128{
129 statCode = input;
130}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNADingfelderChargeIncreaseModel::fParticleChangeForGamma
protected

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