Geant4 11.2.2
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=nullptr, const G4String &nam="DNADingfelderChargeDecreaseModel")
 
 ~G4DNADingfelderChargeDecreaseModel () override=default
 
G4DNADingfelderChargeDecreaseModeloperator= (const G4DNADingfelderChargeDecreaseModel &right)=delete
 
 G4DNADingfelderChargeDecreaseModel (const G4DNADingfelderChargeDecreaseModel &)=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 G4DNADingfelderChargeDecreaseModel.hh.

Constructor & Destructor Documentation

◆ G4DNADingfelderChargeDecreaseModel() [1/2]

G4DNADingfelderChargeDecreaseModel::G4DNADingfelderChargeDecreaseModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "DNADingfelderChargeDecreaseModel" )

Definition at line 46 of file G4DNADingfelderChargeDecreaseModel.cc.

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

◆ ~G4DNADingfelderChargeDecreaseModel()

G4DNADingfelderChargeDecreaseModel::~G4DNADingfelderChargeDecreaseModel ( )
overridedefault

◆ G4DNADingfelderChargeDecreaseModel() [2/2]

G4DNADingfelderChargeDecreaseModel::G4DNADingfelderChargeDecreaseModel ( const G4DNADingfelderChargeDecreaseModel & )
delete

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 208 of file G4DNADingfelderChargeDecreaseModel.cc.

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

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 73 of file G4DNADingfelderChargeDecreaseModel.cc.

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

◆ operator=()

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 279 of file G4DNADingfelderChargeDecreaseModel.cc.

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

◆ SelectStationary()

void G4DNADingfelderChargeDecreaseModel::SelectStationary ( G4bool input)
inline

Definition at line 126 of file G4DNADingfelderChargeDecreaseModel.hh.

127{
128 statCode = input;
129}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNADingfelderChargeDecreaseModel::fParticleChangeForGamma = nullptr
protected

Definition at line 70 of file G4DNADingfelderChargeDecreaseModel.hh.

Referenced by Initialise(), and SampleSecondaries().


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