Geant4 9.6.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)
 
- 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 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 ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., 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 *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
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)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *)
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
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
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void ForceBuildTable (G4bool val)
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

Protected Attributes

G4ParticleChangeForGammafParticleChangeForGamma
 
- Protected Attributes inherited from G4VEmModel
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 

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

Constructor & Destructor Documentation

◆ G4DNADingfelderChargeIncreaseModel()

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

Definition at line 40 of file G4DNADingfelderChargeIncreaseModel.cc.

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

◆ ~G4DNADingfelderChargeIncreaseModel()

G4DNADingfelderChargeIncreaseModel::~G4DNADingfelderChargeIncreaseModel ( )
virtual

Definition at line 67 of file G4DNADingfelderChargeIncreaseModel.cc.

68{}

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

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

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 72 of file G4DNADingfelderChargeIncreaseModel.cc.

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 286 of file G4DNADingfelderChargeIncreaseModel.cc.

291{
292 if (verboseLevel > 3)
293 G4cout << "Calling SampleSecondaries() of G4DNADingfelderChargeIncreaseModel" << G4endl;
294
297
298 G4ParticleDefinition* definition = aDynamicParticle->GetDefinition();
299
300 G4double particleMass = definition->GetPDGMass();
301
302 G4double inK = aDynamicParticle->GetKineticEnergy();
303
304 G4int finalStateIndex = RandomSelect(inK,definition);
305
306 G4int n = NumberOfFinalStates(definition,finalStateIndex);
307
308 G4double outK = inK - IncomingParticleBindingEnergyConstant(definition,finalStateIndex);
309
310 G4DNAGenericIonsManager* instance;
312
313 G4double electronK;
314 if (definition == instance->GetIon("hydrogen")) electronK = inK*electron_mass_c2/proton_mass_c2;
315 else electronK = inK*electron_mass_c2/(particleMass);
316
317 if (outK<0)
318 {
319 G4Exception("G4DNADingfelderChargeIncreaseModel::SampleSecondaries","em0004",
320 FatalException,"Final kinetic energy is negative.");
321 }
322
324
325 (OutgoingParticleDefinition(definition,finalStateIndex),
326 aDynamicParticle->GetMomentumDirection(),
327 outK) ;
328
329 fvect->push_back(dp);
330
331 n = n - 1;
332
333 while (n>0)
334 {
335 n--;
336 fvect->push_back(new G4DynamicParticle
337 (G4Electron::Electron(), aDynamicParticle->GetMomentumDirection(), electronK) );
338 }
339
340}
@ FatalException
@ fStopAndKill
int G4int
Definition: G4Types.hh:66
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition: G4Electron.cc:94
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNADingfelderChargeIncreaseModel::fParticleChangeForGamma
protected

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