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

#include <G4DNADingfelderChargeDecreaseModel.hh>

+ Inheritance diagram for G4DNADingfelderChargeDecreaseModel:

Public Member Functions

 G4DNADingfelderChargeDecreaseModel (const G4ParticleDefinition *p=0, const G4String &nam="DNADingfelderChargeDecreaseModel")
 
virtual ~G4DNADingfelderChargeDecreaseModel ()
 
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 39 of file G4DNADingfelderChargeDecreaseModel.hh.

Constructor & Destructor Documentation

◆ G4DNADingfelderChargeDecreaseModel()

G4DNADingfelderChargeDecreaseModel::G4DNADingfelderChargeDecreaseModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNADingfelderChargeDecreaseModel" 
)

Definition at line 40 of file G4DNADingfelderChargeDecreaseModel.cc.

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

◆ ~G4DNADingfelderChargeDecreaseModel()

G4DNADingfelderChargeDecreaseModel::~G4DNADingfelderChargeDecreaseModel ( )
virtual

Definition at line 70 of file G4DNADingfelderChargeDecreaseModel.cc.

71{
72}

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 209 of file G4DNADingfelderChargeDecreaseModel.cc.

214{
215 if (verboseLevel > 3)
216 {
217 G4cout
218 << "Calling CrossSectionPerVolume() of G4DNADingfelderChargeDecreaseModel"
219 << G4endl;
220 }
221
222 // Calculate total cross section for model
223
224 G4DNAGenericIonsManager *instance;
226
227 if (
228 particleDefinition != G4Proton::ProtonDefinition()
229 &&
230 particleDefinition != instance->GetIon("alpha++")
231 &&
232 particleDefinition != instance->GetIon("alpha+")
233 )
234
235 return 0;
236
237 G4double lowLim = 0;
238 G4double highLim = 0;
239 G4double crossSection = 0.;
240
241 G4double waterDensity = (*fpMolWaterDensity)[material->GetIndex()];
242
243 const G4String& particleName = particleDefinition->GetParticleName();
244
245 std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
246 pos1 = lowEnergyLimit.find(particleName);
247
248 if (pos1 != lowEnergyLimit.end())
249 {
250 lowLim = pos1->second;
251 }
252
253 std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
254 pos2 = highEnergyLimit.find(particleName);
255
256 if (pos2 != highEnergyLimit.end())
257 {
258 highLim = pos2->second;
259 }
260
261 if (k >= lowLim && k <= highLim)
262 {
263 crossSection = Sum(k,particleDefinition);
264 }
265
266 if (verboseLevel > 2)
267 {
268 G4cout << "_______________________________________" << G4endl;
269 G4cout << "G4DNADingfelderChargeDecreaeModel" << G4endl;
270 G4cout << "Kinetic energy(eV)=" << k/eV << "particle :" << particleName << G4endl;
271 G4cout << "Cross section per water molecule (cm^2)=" << crossSection/cm/cm << G4endl;
272 G4cout << "Cross section per water molecule (cm^-1)=" << crossSection*
273 waterDensity/(1./cm) << G4endl;
274 // material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
275 }
276
277 return crossSection*waterDensity;
278 //return crossSection*material->GetAtomicNumDensityVector()[1];
279
280}
double G4double
Definition: G4Types.hh:83
static G4DNAGenericIonsManager * Instance(void)
G4ParticleDefinition * GetIon(const G4String &name)
size_t GetIndex() const
Definition: G4Material.hh:258
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:87

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 76 of file G4DNADingfelderChargeDecreaseModel.cc.

78{
79
80 if (verboseLevel > 3)
81 {
82 G4cout << "Calling G4DNADingfelderChargeDecreaseModel::Initialise()"
83 << G4endl;
84 }
85
86 // Energy limits
87
91 G4ParticleDefinition* alphaPlusPlusDef = instance->GetIon("alpha++");
92 G4ParticleDefinition* alphaPlusDef = instance->GetIon("alpha+");
93
95 G4String alphaPlusPlus;
96 G4String alphaPlus;
97
98 // Limits
99
100 proton = protonDef->GetParticleName();
101 lowEnergyLimit[proton] = 100. * eV;
102 highEnergyLimit[proton] = 100. * MeV;
103
104 alphaPlusPlus = alphaPlusPlusDef->GetParticleName();
105 lowEnergyLimit[alphaPlusPlus] = 1. * keV;
106 highEnergyLimit[alphaPlusPlus] = 400. * MeV;
107
108 alphaPlus = alphaPlusDef->GetParticleName();
109 lowEnergyLimit[alphaPlus] = 1. * keV;
110 highEnergyLimit[alphaPlus] = 400. * MeV;
111
112 //
113
114 if (particle==protonDef)
115 {
116 SetLowEnergyLimit(lowEnergyLimit[proton]);
117 SetHighEnergyLimit(highEnergyLimit[proton]);
118 }
119
120 if (particle==alphaPlusPlusDef)
121 {
122 SetLowEnergyLimit(lowEnergyLimit[alphaPlusPlus]);
123 SetHighEnergyLimit(highEnergyLimit[alphaPlusPlus]);
124 }
125
126 if (particle==alphaPlusDef)
127 {
128 SetLowEnergyLimit(lowEnergyLimit[alphaPlus]);
129 SetHighEnergyLimit(highEnergyLimit[alphaPlus]);
130 }
131
132 // Final state
133
134 //PROTON
135 f0[0][0]=1.;
136 a0[0][0]=-0.180;
137 a1[0][0]=-3.600;
138 b0[0][0]=-18.22;
139 b1[0][0]=-1.997;
140 c0[0][0]=0.215;
141 d0[0][0]=3.550;
142 x0[0][0]=3.450;
143 x1[0][0]=5.251;
144
145 numberOfPartialCrossSections[0] = 1;
146
147 //ALPHA++
148 f0[0][1]=1.; a0[0][1]=0.95;
149 a1[0][1]=-2.75;
150 b0[0][1]=-23.00;
151 c0[0][1]=0.215;
152 d0[0][1]=2.95;
153 x0[0][1]=3.50;
154
155 f0[1][1]=1.;
156 a0[1][1]=0.95;
157 a1[1][1]=-2.75;
158 b0[1][1]=-23.73;
159 c0[1][1]=0.250;
160 d0[1][1]=3.55;
161 x0[1][1]=3.72;
162
163 x1[0][1]=-1.;
164 b1[0][1]=-1.;
165
166 x1[1][1]=-1.;
167 b1[1][1]=-1.;
168
169 numberOfPartialCrossSections[1] = 2;
170
171 // ALPHA+
172 f0[0][2]=1.;
173 a0[0][2]=0.65;
174 a1[0][2]=-2.75;
175 b0[0][2]=-21.81;
176 c0[0][2]=0.232;
177 d0[0][2]=2.95;
178 x0[0][2]=3.53;
179
180 x1[0][2]=-1.;
181 b1[0][2]=-1.;
182
183 numberOfPartialCrossSections[2] = 1;
184
185 //
186
187 if( verboseLevel>0 )
188 {
189 G4cout << "Dingfelder charge decrease model is initialized " << G4endl
190 << "Energy range: "
191 << LowEnergyLimit() / keV << " keV - "
192 << HighEnergyLimit() / MeV << " MeV for "
193 << particle->GetParticleName()
194 << G4endl;
195 }
196
197 // Initialize water density pointer
199
200 if (isInitialised)
201 { return;}
203 isInitialised = true;
204
205}
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 G4DNADingfelderChargeDecreaseModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple ,
const G4DynamicParticle aDynamicParticle,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 284 of file G4DNADingfelderChargeDecreaseModel.cc.

290{
291 if (verboseLevel > 3)
292 {
293 G4cout
294 << "Calling SampleSecondaries() of G4DNADingfelderChargeDecreaseModel"
295 << G4endl;
296 }
297
298 G4double inK = aDynamicParticle->GetKineticEnergy();
299
300 G4ParticleDefinition* definition = aDynamicParticle->GetDefinition();
301
302 G4double particleMass = definition->GetPDGMass();
303
304 G4int finalStateIndex = RandomSelect(inK,definition);
305
306 G4int n = NumberOfFinalStates(definition, finalStateIndex);
307 G4double waterBindingEnergy = WaterBindingEnergyConstant(definition, finalStateIndex);
308 G4double outgoingParticleBindingEnergy = OutgoingParticleBindingEnergyConstant(definition, finalStateIndex);
309
310 G4double outK = 0.;
311
312 if (definition==G4Proton::Proton())
313 {
314 if (!statCode) outK = inK - n*(inK*electron_mass_c2/proton_mass_c2)
315 - waterBindingEnergy + outgoingParticleBindingEnergy;
316 else outK = inK;
317 }
318
319 else
320 {
321 if (!statCode) outK = inK - n*(inK*electron_mass_c2/particleMass)
322 - waterBindingEnergy + outgoingParticleBindingEnergy;
323 else outK = inK;
324 }
325
326 if (outK<0)
327 {
328 G4Exception("G4DNADingfelderChargeDecreaseModel::SampleSecondaries","em0004",
329 FatalException,"Final kinetic energy is negative.");
330 }
331
333
334 if (!statCode) fParticleChangeForGamma->ProposeLocalEnergyDeposit(waterBindingEnergy);
335
336 else
337
338 {
339 if (definition==G4Proton::Proton())
340 fParticleChangeForGamma->ProposeLocalEnergyDeposit(n*(inK*electron_mass_c2/proton_mass_c2)
341 + waterBindingEnergy - outgoingParticleBindingEnergy);
342 else
343 fParticleChangeForGamma->ProposeLocalEnergyDeposit(n*(inK*electron_mass_c2/particleMass)
344 + waterBindingEnergy - outgoingParticleBindingEnergy);
345 }
346
347 G4DynamicParticle* dp = new G4DynamicParticle (OutgoingParticleDefinition(definition, finalStateIndex),
348 aDynamicParticle->GetMomentumDirection(),
349 outK);
350 fvect->push_back(dp);
351
352 const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
354 1,
355 theIncomingTrack);
356
357}
@ eIonizedMolecule
@ 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
static G4DNAChemistryManager * Instance()
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
const G4Track * GetCurrentTrack() const
static G4Proton * Proton()
Definition: G4Proton.cc:92
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

◆ SelectStationary()

void G4DNADingfelderChargeDecreaseModel::SelectStationary ( G4bool  input)
inline

Definition at line 121 of file G4DNADingfelderChargeDecreaseModel.hh.

122{
123 statCode = input;
124}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNADingfelderChargeDecreaseModel::fParticleChangeForGamma
protected

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