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

#include <G4DNARuddIonisationExtendedModel.hh>

+ Inheritance diagram for G4DNARuddIonisationExtendedModel:

Public Member Functions

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

Constructor & Destructor Documentation

◆ G4DNARuddIonisationExtendedModel()

G4DNARuddIonisationExtendedModel::G4DNARuddIonisationExtendedModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNARuddIonisationExtendedModel" 
)

Definition at line 49 of file G4DNARuddIonisationExtendedModel.cc.

51:G4VEmModel(nam),isInitialised(false)
52{
53 // nistwater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER");
54 fpWaterDensity = 0;
55
56 slaterEffectiveCharge[0]=0.;
57 slaterEffectiveCharge[1]=0.;
58 slaterEffectiveCharge[2]=0.;
59 sCoefficient[0]=0.;
60 sCoefficient[1]=0.;
61 sCoefficient[2]=0.;
62
63 lowEnergyLimitForA[1] = 0 * eV;
64 lowEnergyLimitForA[2] = 0 * eV;
65 lowEnergyLimitForA[3] = 0 * eV;
66 lowEnergyLimitOfModelForA[1] = 100 * eV;
67 lowEnergyLimitOfModelForA[4] = 1 * keV;
68 lowEnergyLimitOfModelForA[5] = 0.5 * MeV; // For A = 3 or above, limit is MeV/uma
69 killBelowEnergyForA[1] = lowEnergyLimitOfModelForA[1];
70 killBelowEnergyForA[4] = lowEnergyLimitOfModelForA[4];
71 killBelowEnergyForA[5] = lowEnergyLimitOfModelForA[5];
72
73 verboseLevel= 0;
74 // Verbosity scale:
75 // 0 = nothing
76 // 1 = warning for energy non-conservation
77 // 2 = details of energy budget
78 // 3 = calculation of cross sections, file openings, sampling of atoms
79 // 4 = entering in methods
80
81 if( verboseLevel>0 )
82 {
83 G4cout << "Rudd ionisation model is constructed " << G4endl;
84 }
85
86 // Define default angular generator
88
89 // Mark this model as "applicable" for atomic deexcitation
91 fAtomDeexcitation = 0;
93
94 // Selection of stationary mode
95
96 statCode = false;
97}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:813
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:618

◆ ~G4DNARuddIonisationExtendedModel()

G4DNARuddIonisationExtendedModel::~G4DNARuddIonisationExtendedModel ( )
virtual

Definition at line 101 of file G4DNARuddIonisationExtendedModel.cc.

102{
103 // Cross section
104
105 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
106 for (pos = tableData.begin(); pos != tableData.end(); ++pos)
107 {
108 G4DNACrossSectionDataSet* table = pos->second;
109 delete table;
110 }
111
112 // The following removal is forbidden G4VEnergyLossModel takes care of deletion
113 // however coverity will signal this as an error
114 // if (fAtomDeexcitation) {delete fAtomDeexcitation;}
115
116}

Member Function Documentation

◆ CrossSectionPerVolume()

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

Reimplemented from G4VEmModel.

Definition at line 535 of file G4DNARuddIonisationExtendedModel.cc.

540{
541 //SI: particleDefinition->GetParticleName() is for eg. Fe56
542 // particleDefinition->GetPDGMass() is correct
543 // particleDefinition->GetAtomicNumber() is correct
544
545 if (verboseLevel > 3)
546 G4cout << "Calling CrossSectionPerVolume() of G4DNARuddIonisationExtendedModel" << G4endl;
547
548 // Calculate total cross section for model
549
550 G4DNAGenericIonsManager *instance;
552
553 if (
554 particleDefinition != G4Proton::ProtonDefinition()
555 &&
556 particleDefinition != instance->GetIon("hydrogen")
557 &&
558 particleDefinition != instance->GetIon("alpha++")
559 &&
560 particleDefinition != instance->GetIon("alpha+")
561 &&
562 particleDefinition != instance->GetIon("helium")
563 &&
564 // SI
565 //particleDefinition != instance->GetIon("carbon")
566 //&&
567 //particleDefinition != instance->GetIon("nitrogen")
568 //&&
569 //particleDefinition != instance->GetIon("oxygen")
570 //&&
571 //particleDefinition != instance->GetIon("iron")
572 particleDefinition != G4IonTable::GetIonTable()->GetIon(3,7)
573 &&
574 particleDefinition != G4IonTable::GetIonTable()->GetIon(4,9)
575 &&
576 particleDefinition != G4IonTable::GetIonTable()->GetIon(5,11)
577 &&
578 particleDefinition != G4IonTable::GetIonTable()->GetIon(6,12)
579 &&
580 particleDefinition != G4IonTable::GetIonTable()->GetIon(7,14)
581 &&
582 particleDefinition != G4IonTable::GetIonTable()->GetIon(8,16)
583 &&
584 particleDefinition != G4IonTable::GetIonTable()->GetIon(14,28)
585 &&
586 particleDefinition != G4IonTable::GetIonTable()->GetIon(26,56)
587 //
588 )
589
590 return 0;
591
592 G4double lowLim = 0;
593
594 if ( particleDefinition == G4Proton::ProtonDefinition()
595 || particleDefinition == instance->GetIon("hydrogen")
596 )
597
598 lowLim = lowEnergyLimitOfModelForA[1];
599
600 else if ( particleDefinition == instance->GetIon("alpha++")
601 || particleDefinition == instance->GetIon("alpha+")
602 || particleDefinition == instance->GetIon("helium")
603 )
604
605 lowLim = lowEnergyLimitOfModelForA[4];
606
607 else lowLim = lowEnergyLimitOfModelForA[5];
608
609 G4double highLim = 0;
610 G4double sigma=0;
611
612
613 G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
614
615 const G4String& particleName = particleDefinition->GetParticleName();
616
617 std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
618 pos2 = highEnergyLimit.find(particleName);
619
620 if (pos2 != highEnergyLimit.end())
621 {
622 highLim = pos2->second;
623 }
624
625 if (k <= highLim)
626 {
627
628 //SI : XS must not be zero otherwise sampling of secondaries method ignored
629
630 if (k < lowLim) k = lowLim;
631
632 //
633
634 std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
635 pos = tableData.find(particleName);
636
637 if (pos != tableData.end())
638 {
639 G4DNACrossSectionDataSet* table = pos->second;
640 if (table != 0)
641 {
642 sigma = table->FindValue(k);
643 }
644 }
645 else
646 {
647 G4Exception("G4DNARuddIonisationExtendedModel::CrossSectionPerVolume","em0002",
648 FatalException,"Model not applicable to particle type.");
649 }
650
651 } // if (k >= lowLim && k < highLim)
652
653 if (verboseLevel > 2)
654 {
655 G4cout << "__________________________________" << G4endl;
656 G4cout << "G4DNARuddIonisationExtendedModel - XS INFO START" << G4endl;
657 G4cout << "Kinetic energy(eV)=" << k/eV << " particle : " << particleDefinition->GetParticleName() << G4endl;
658 G4cout << "Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
659 G4cout << "Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
660 //G4cout << " - Cross section per water molecule (cm^-1)="
661 //<< sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
662 G4cout << "G4DNARuddIonisationExtendedModel - XS INFO END" << G4endl;
663
664 }
665
666 return sigma*waterDensity;
667
668}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
double G4double
Definition: G4Types.hh:83
virtual G4double FindValue(G4double e, G4int componentId=0) const
static G4DNAGenericIonsManager * Instance(void)
G4ParticleDefinition * GetIon(const G4String &name)
G4ParticleDefinition * GetIon(G4int Z, G4int A, G4int lvl=0)
Definition: G4IonTable.cc:522
static G4IonTable * GetIonTable()
Definition: G4IonTable.cc:170
size_t GetIndex() const
Definition: G4Material.hh:258
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:87

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 120 of file G4DNARuddIonisationExtendedModel.cc.

122{
123 if (verboseLevel > 3)
124 G4cout << "Calling G4DNARuddIonisationExtendedModel::Initialise()" << G4endl;
125
126 // Energy limits
127
128 G4String fileProton("dna/sigma_ionisation_p_rudd");
129 G4String fileHydrogen("dna/sigma_ionisation_h_rudd");
130 G4String fileAlphaPlusPlus("dna/sigma_ionisation_alphaplusplus_rudd");
131 G4String fileAlphaPlus("dna/sigma_ionisation_alphaplus_rudd");
132 G4String fileHelium("dna/sigma_ionisation_he_rudd");
133 G4String fileLithium("dna/sigma_ionisation_li_rudd");
134 G4String fileBeryllium("dna/sigma_ionisation_be_rudd");
135 G4String fileBoron("dna/sigma_ionisation_b_rudd");
136 G4String fileCarbon("dna/sigma_ionisation_c_rudd");
137 G4String fileNitrogen("dna/sigma_ionisation_n_rudd");
138 G4String fileOxygen("dna/sigma_ionisation_o_rudd");
139 G4String fileSilicon("dna/sigma_ionisation_si_rudd");
140 G4String fileIron("dna/sigma_ionisation_fe_rudd");
141
142 G4DNAGenericIonsManager *instance;
145 G4ParticleDefinition* hydrogenDef = instance->GetIon("hydrogen");
146 G4ParticleDefinition* alphaPlusPlusDef = instance->GetIon("alpha++");
147 G4ParticleDefinition* alphaPlusDef = instance->GetIon("alpha+");
148 G4ParticleDefinition* heliumDef = instance->GetIon("helium");
149
150 //G4ParticleDefinition* carbonDef = instance->GetIon("carbon");
151 //G4ParticleDefinition* nitrogenDef = instance->GetIon("nitrogen");
152 //G4ParticleDefinition* oxygenDef = instance->GetIon("oxygen");
153 //G4ParticleDefinition* siliconDef = instance->GetIon("silicon");
154 //G4ParticleDefinition* ironDef = instance->GetIon("iron");
163 //
164
166 G4String hydrogen;
167 G4String alphaPlusPlus;
168 G4String alphaPlus;
169 G4String helium;
170 G4String lithium;
171 G4String beryllium;
172 G4String boron;
173 G4String carbon;
174 G4String nitrogen;
175 G4String oxygen;
176 G4String silicon;
177 G4String iron;
178
179 G4double scaleFactor = 1 * m*m;
180
181 // LIMITS AND DATA
182
183 // **********************************************************************************************
184
185 proton = protonDef->GetParticleName();
186 tableFile[proton] = fileProton;
187 lowEnergyLimit[proton] = lowEnergyLimitForA[1];
188 highEnergyLimit[proton] = 500. * keV;
189
190 // Cross section
191
193 eV,
194 scaleFactor );
195 tableProton->LoadData(fileProton);
196 tableData[proton] = tableProton;
197
198 // **********************************************************************************************
199
200 hydrogen = hydrogenDef->GetParticleName();
201 tableFile[hydrogen] = fileHydrogen;
202
203 lowEnergyLimit[hydrogen] = lowEnergyLimitForA[1];
204 highEnergyLimit[hydrogen] = 100. * MeV;
205
206 // Cross section
207
209 eV,
210 scaleFactor );
211 tableHydrogen->LoadData(fileHydrogen);
212
213 tableData[hydrogen] = tableHydrogen;
214
215 // **********************************************************************************************
216
217 alphaPlusPlus = alphaPlusPlusDef->GetParticleName();
218 tableFile[alphaPlusPlus] = fileAlphaPlusPlus;
219
220 lowEnergyLimit[alphaPlusPlus] = lowEnergyLimitForA[4];
221 highEnergyLimit[alphaPlusPlus] = 400. * MeV;
222
223 // Cross section
224
226 eV,
227 scaleFactor );
228 tableAlphaPlusPlus->LoadData(fileAlphaPlusPlus);
229
230 tableData[alphaPlusPlus] = tableAlphaPlusPlus;
231
232 // **********************************************************************************************
233
234 alphaPlus = alphaPlusDef->GetParticleName();
235 tableFile[alphaPlus] = fileAlphaPlus;
236
237 lowEnergyLimit[alphaPlus] = lowEnergyLimitForA[4];
238 highEnergyLimit[alphaPlus] = 400. * MeV;
239
240 // Cross section
241
243 eV,
244 scaleFactor );
245 tableAlphaPlus->LoadData(fileAlphaPlus);
246 tableData[alphaPlus] = tableAlphaPlus;
247
248 // **********************************************************************************************
249
250 helium = heliumDef->GetParticleName();
251 tableFile[helium] = fileHelium;
252
253 lowEnergyLimit[helium] = lowEnergyLimitForA[4];
254 highEnergyLimit[helium] = 400. * MeV;
255
256 // Cross section
257
259 eV,
260 scaleFactor );
261 tableHelium->LoadData(fileHelium);
262 tableData[helium] = tableHelium;
263
264 // **********************************************************************************************
265
266 lithium = lithiumDef->GetParticleName();
267 tableFile[lithium] = fileLithium;
268
269 //SI
270 //lowEnergyLimit[carbon] = lowEnergyLimitForA[5] * particle->GetAtomicMass();
271 //highEnergyLimit[carbon] = 1e6* particle->GetAtomicMass() * MeV;
272 lowEnergyLimit[lithium] = 0.5*7*MeV;
273 highEnergyLimit[lithium] = 1e6*7*MeV;
274 //
275
276 // Cross section
277
279 eV,
280 scaleFactor );
281 tableLithium->LoadData(fileLithium);
282 tableData[lithium] = tableLithium;
283
284 // **********************************************************************************************
285
286 beryllium = berylliumDef->GetParticleName();
287 tableFile[beryllium] = fileBeryllium;
288
289 //SI
290 //lowEnergyLimit[carbon] = lowEnergyLimitForA[5] * particle->GetAtomicMass();
291 //highEnergyLimit[carbon] = 1e6* particle->GetAtomicMass() * MeV;
292 lowEnergyLimit[beryllium] = 0.5*9*MeV;
293 highEnergyLimit[beryllium] = 1e6*9*MeV;
294 //
295
296 // Cross section
297
299 eV,
300 scaleFactor );
301 tableBeryllium->LoadData(fileBeryllium);
302 tableData[beryllium] = tableBeryllium;
303
304 // **********************************************************************************************
305
306 boron = boronDef->GetParticleName();
307 tableFile[boron] = fileBoron;
308
309 //SI
310 //lowEnergyLimit[carbon] = lowEnergyLimitForA[5] * particle->GetAtomicMass();
311 //highEnergyLimit[carbon] = 1e6* particle->GetAtomicMass() * MeV;
312 lowEnergyLimit[boron] = 0.5*11*MeV;
313 highEnergyLimit[boron] = 1e6*11*MeV;
314 //
315
316 // Cross section
317
319 eV,
320 scaleFactor );
321 tableBoron->LoadData(fileBoron);
322 tableData[boron] = tableBoron;
323
324 // **********************************************************************************************
325
326 carbon = carbonDef->GetParticleName();
327 tableFile[carbon] = fileCarbon;
328
329 //SI
330 //lowEnergyLimit[carbon] = lowEnergyLimitForA[5] * particle->GetAtomicMass();
331 //highEnergyLimit[carbon] = 1e6* particle->GetAtomicMass() * MeV;
332 lowEnergyLimit[carbon] = 0.5*12*MeV;
333 highEnergyLimit[carbon] = 1e6*12*MeV;
334 //
335
336 // Cross section
337
339 eV,
340 scaleFactor );
341 tableCarbon->LoadData(fileCarbon);
342 tableData[carbon] = tableCarbon;
343
344 // **********************************************************************************************
345
346 oxygen = oxygenDef->GetParticleName();
347 tableFile[oxygen] = fileOxygen;
348
349 //SI
350 //lowEnergyLimit[oxygen] = lowEnergyLimitForA[5]* particle->GetAtomicMass();
351 //highEnergyLimit[oxygen] = 1e6* particle->GetAtomicMass()* MeV;
352 lowEnergyLimit[oxygen] = 0.5*16*MeV;
353 highEnergyLimit[oxygen] = 1e6*16*MeV;
354 //
355
356 // Cross section
357
359 eV,
360 scaleFactor );
361 tableOxygen->LoadData(fileOxygen);
362 tableData[oxygen] = tableOxygen;
363
364 // **********************************************************************************************
365
366 nitrogen = nitrogenDef->GetParticleName();
367 tableFile[nitrogen] = fileNitrogen;
368
369 //SI
370 //lowEnergyLimit[nitrogen] = lowEnergyLimitForA[5]* particle->GetAtomicMass();
371 //highEnergyLimit[nitrogen] = 1e6* particle->GetAtomicMass()* MeV;
372 lowEnergyLimit[nitrogen] = 0.5*14*MeV;
373 highEnergyLimit[nitrogen] = 1e6*14*MeV;
374 //
375
376 // Cross section
377
379 eV,
380 scaleFactor );
381 tableNitrogen->LoadData(fileNitrogen);
382 tableData[nitrogen] = tableNitrogen;
383
384 // **********************************************************************************************
385
386 silicon = siliconDef->GetParticleName();
387 tableFile[silicon] = fileSilicon;
388
389 //lowEnergyLimit[silicon] = lowEnergyLimitForA[5]* particle->GetAtomicMass();
390 //highEnergyLimit[silicon] = 1e6* particle->GetAtomicMass()* MeV;
391 lowEnergyLimit[silicon] = 0.5*28*MeV;
392 highEnergyLimit[silicon] = 1e6*28*MeV;
393 //
394
395 // Cross section
396
398 eV,
399 scaleFactor );
400 tableSilicon->LoadData(fileSilicon);
401 tableData[silicon] = tableSilicon;
402
403 // **********************************************************************************************
404
405 iron = ironDef->GetParticleName();
406 tableFile[iron] = fileIron;
407
408 //SI
409 //lowEnergyLimit[iron] = lowEnergyLimitForA[5]* particle->GetAtomicMass();
410 //highEnergyLimit[iron] = 1e6* particle->GetAtomicMass()* MeV;
411 lowEnergyLimit[iron] = 0.5*56*MeV;
412 highEnergyLimit[iron] = 1e6*56*MeV;
413 //
414
415 // Cross section
416
418 eV,
419 scaleFactor );
420 tableIron->LoadData(fileIron);
421 tableData[iron] = tableIron;
422
423 // **********************************************************************************************
424
425 // SI: not anymore
426 // ZF Following lines can be replaced by:
427 // SetLowEnergyLimit(lowEnergyLimit[particle->GetParticleName()]);
428 // SetHighEnergyLimit(highEnergyLimit[particle->GetParticleName()]);
429 // at least for HZE
430
431 if (particle==protonDef)
432 {
433 SetLowEnergyLimit(lowEnergyLimit[proton]);
434 SetHighEnergyLimit(highEnergyLimit[proton]);
435 }
436
437 if (particle==hydrogenDef)
438 {
439 SetLowEnergyLimit(lowEnergyLimit[hydrogen]);
440 SetHighEnergyLimit(highEnergyLimit[hydrogen]);
441 }
442
443 if (particle==heliumDef)
444 {
445 SetLowEnergyLimit(lowEnergyLimit[helium]);
446 SetHighEnergyLimit(highEnergyLimit[helium]);
447 }
448
449 if (particle==alphaPlusDef)
450 {
451 SetLowEnergyLimit(lowEnergyLimit[alphaPlus]);
452 SetHighEnergyLimit(highEnergyLimit[alphaPlus]);
453 }
454
455 if (particle==alphaPlusPlusDef)
456 {
457 SetLowEnergyLimit(lowEnergyLimit[alphaPlusPlus]);
458 SetHighEnergyLimit(highEnergyLimit[alphaPlusPlus]);
459 }
460
461 if (particle==lithiumDef)
462 {
463 SetLowEnergyLimit(lowEnergyLimit[lithium]);
464 SetHighEnergyLimit(highEnergyLimit[lithium]);
465 }
466
467 if (particle==berylliumDef)
468 {
469 SetLowEnergyLimit(lowEnergyLimit[beryllium]);
470 SetHighEnergyLimit(highEnergyLimit[beryllium]);
471 }
472
473 if (particle==boronDef)
474 {
475 SetLowEnergyLimit(lowEnergyLimit[boron]);
476 SetHighEnergyLimit(highEnergyLimit[boron]);
477 }
478
479 if (particle==carbonDef)
480 {
481 SetLowEnergyLimit(lowEnergyLimit[carbon]);
482 SetHighEnergyLimit(highEnergyLimit[carbon]);
483 }
484
485 if (particle==nitrogenDef)
486 {
487 SetLowEnergyLimit(lowEnergyLimit[nitrogen]);
488 SetHighEnergyLimit(highEnergyLimit[nitrogen]);
489 }
490
491 if (particle==oxygenDef)
492 {
493 SetLowEnergyLimit(lowEnergyLimit[oxygen]);
494 SetHighEnergyLimit(highEnergyLimit[oxygen]);
495 }
496
497 if (particle==siliconDef)
498 {
499 SetLowEnergyLimit(lowEnergyLimit[silicon]);
500 SetHighEnergyLimit(highEnergyLimit[silicon]);
501 }
502
503 if (particle==ironDef)
504 {
505 SetLowEnergyLimit(lowEnergyLimit[iron]);
506 SetHighEnergyLimit(highEnergyLimit[iron]);
507 }
508
509 //----------------------------------------------------------------------
510
511 if( verboseLevel>0 )
512 {
513 G4cout << "Rudd ionisation model is initialized " << G4endl
514 << "Energy range: "
515 << LowEnergyLimit() / eV << " eV - "
516 << HighEnergyLimit() / keV << " keV for "
517 << particle->GetParticleName()
518 << G4endl;
519 }
520
521 // Initialize water density pointer
523
524 //
525
526 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
527
528 if (isInitialised) { return; }
530 isInitialised = true;
531}
virtual G4bool LoadData(const G4String &argFileName)
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 G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
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 G4DNARuddIonisationExtendedModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle particle,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 672 of file G4DNARuddIonisationExtendedModel.cc.

677{
678 //SI: particle->GetDefinition()->GetParticleName() is for eg. Fe56
679 // particle->GetDefinition()->GetPDGMass() is correct
680 // particle->GetDefinition()->GetAtomicNumber() is correct
681 // particle->GetDefinition()->GetAtomicMass() is correct
682
683 if (verboseLevel > 3)
684 G4cout << "Calling SampleSecondaries() of G4DNARuddIonisationExtendedModel" << G4endl;
685
686 G4double lowLim = 0;
687 G4double highLim = 0;
688
689 // ZF: the following line summarizes the commented part
690
691 if(particle->GetDefinition()->GetAtomicMass() <= 4) lowLim = killBelowEnergyForA[particle->GetDefinition()->GetAtomicMass()];
692
693 else lowLim = killBelowEnergyForA[5]*particle->GetDefinition()->GetAtomicMass();
694
695 /*
696
697 if(particle->GetDefinition()->GetAtomicMass() >= 5) lowLim = killBelowEnergyForA[5]*particle->GetDefinition()->GetAtomicMass();
698
699 if ( particle->GetDefinition() == G4Proton::ProtonDefinition()
700 || particle->GetDefinition() == instance->GetIon("hydrogen")
701 )
702
703 lowLim = killBelowEnergyForA[1];
704
705 if ( particle->GetDefinition() == instance->GetIon("alpha++")
706 || particle->GetDefinition() == instance->GetIon("alpha+")
707 || particle->GetDefinition() == instance->GetIon("helium")
708 )
709
710 lowLim = killBelowEnergyForA[4];
711
712 */
713 //
714
715 G4double k = particle->GetKineticEnergy();
716
717 const G4String& particleName = particle->GetDefinition()->GetParticleName();
718
719 // SI - the following is useless since lowLim is already defined
720 /*
721 std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
722 pos1 = lowEnergyLimit.find(particleName);
723
724 if (pos1 != lowEnergyLimit.end())
725 {
726 lowLim = pos1->second;
727 }
728 */
729
730 std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
731 pos2 = highEnergyLimit.find(particleName);
732
733 if (pos2 != highEnergyLimit.end()) highLim = pos2->second;
734
735 if (k >= lowLim && k <= highLim)
736
737 // SI: no strict limits, like in the non extended version of the model
738 {
739 G4ParticleDefinition* definition = particle->GetDefinition();
740 G4ParticleMomentum primaryDirection = particle->GetMomentumDirection();
741 /*
742 G4double particleMass = definition->GetPDGMass();
743 G4double totalEnergy = k + particleMass;
744 G4double pSquare = k*(totalEnergy+particleMass);
745 G4double totalMomentum = std::sqrt(pSquare);
746 */
747
748 G4int ionizationShell = RandomSelect(k,particleName);
749
750 // sample deexcitation
751 // here we assume that H_{2}O electronic levels are the same as Oxygen.
752 // this can be considered true with a rough 10% error in energy on K-shell,
753
755 bindingEnergy = waterStructure.IonisationEnergy(ionizationShell);
756
757 //SI: additional protection if tcs interpolation method is modified
758 if (k<bindingEnergy) return;
759 //
760
761 G4double secondaryKinetic = RandomizeEjectedElectronEnergy(definition,k,ionizationShell);
762
763 G4int Z = 8;
764
765 G4ThreeVector deltaDirection =
766 GetAngularDistribution()->SampleDirectionForShell(particle, secondaryKinetic,
767 Z, ionizationShell,
768 couple->GetMaterial());
769
770 G4DynamicParticle* dp = new G4DynamicParticle (G4Electron::Electron(),deltaDirection,secondaryKinetic) ;
771 fvect->push_back(dp);
772
774
775 // SI: the following lines are not needed anymore
776 /*
777 G4double cosTheta = 0.;
778 G4double phi = 0.;
779 RandomizeEjectedElectronDirection(definition, k,secondaryKinetic, cosTheta, phi, ionizationShell);
780
781 G4double sinTheta = std::sqrt(1.-cosTheta*cosTheta);
782 G4double dirX = sinTheta*std::cos(phi);
783 G4double dirY = sinTheta*std::sin(phi);
784 G4double dirZ = cosTheta;
785 G4ThreeVector deltaDirection(dirX,dirY,dirZ);
786 deltaDirection.rotateUz(primaryDirection);
787 */
788
789 // Ignored for ions on electrons
790 /*
791 G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 ));
792
793 G4double finalPx = totalMomentum*primaryDirection.x() - deltaTotalMomentum*deltaDirection.x();
794 G4double finalPy = totalMomentum*primaryDirection.y() - deltaTotalMomentum*deltaDirection.y();
795 G4double finalPz = totalMomentum*primaryDirection.z() - deltaTotalMomentum*deltaDirection.z();
796 G4double finalMomentum = std::sqrt(finalPx*finalPx+finalPy*finalPy+finalPz*finalPz);
797 finalPx /= finalMomentum;
798 finalPy /= finalMomentum;
799 finalPz /= finalMomentum;
800
801 G4ThreeVector direction;
802 direction.set(finalPx,finalPy,finalPz);
803
804 fParticleChangeForGamma->ProposeMomentumDirection(direction.unit()) ;
805 */
806
807 size_t secNumberInit = 0;// need to know at a certain point the energy of secondaries
808 size_t secNumberFinal = 0;// So I'll make the diference and then sum the energies
809
810 G4double scatteredEnergy = k-bindingEnergy-secondaryKinetic;
811
812 // SI: only atomic deexcitation from K shell is considered
813 if(fAtomDeexcitation && ionizationShell == 4)
814 {
815 const G4AtomicShell* shell
816 = fAtomDeexcitation->GetAtomicShell(Z, G4AtomicShellEnumerator(0));
817 secNumberInit = fvect->size();
818 fAtomDeexcitation->GenerateParticles(fvect, shell, Z, 0, 0);
819 secNumberFinal = fvect->size();
820
821 if(secNumberFinal > secNumberInit)
822 {
823 for (size_t i=secNumberInit; i<secNumberFinal; ++i)
824 {
825 //Check if there is enough residual energy
826 if (bindingEnergy >= ((*fvect)[i])->GetKineticEnergy())
827 {
828 //Ok, this is a valid secondary: keep it
829 bindingEnergy -= ((*fvect)[i])->GetKineticEnergy();
830 }
831 else
832 {
833 //Invalid secondary: not enough energy to create it!
834 //Keep its energy in the local deposit
835 delete (*fvect)[i];
836 (*fvect)[i]=0;
837 }
838 }
839 }
840
841 }
842
843 //This should never happen
844 if(bindingEnergy < 0.0)
845 G4Exception("G4DNAEmfietzoglouIonisatioModel1::SampleSecondaries()",
846 "em2050",FatalException,"Negative local energy deposit");
847
848 //bindingEnergy has been decreased
849 //by the amount of energy taken away by deexc. products
850 if (!statCode)
851 {
854 }
855 else
856 {
859 }
860
861 // TEST //////////////////////////
862 // if (secondaryKinetic<0) abort();
863 // if (scatteredEnergy<0) abort();
864 // if (k-scatteredEnergy-secondaryKinetic-deexSecEnergy<0) abort();
865 // if (k-scatteredEnergy<0) abort();
866 /////////////////////////////////
867
868 const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
870 ionizationShell,
871 theIncomingTrack);
872 }
873
874 // SI - not useful since low energy of model is 0 eV
875
876 if (k < lowLim)
877 {
881 }
882
883}
G4AtomicShellEnumerator
@ eIonizedMolecule
@ 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 G4Electron * Electron()
Definition: G4Electron.cc:93
const G4Material * GetMaterial() const
const G4Track * GetCurrentTrack() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4int GetAtomicMass() const
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
virtual G4ThreeVector & SampleDirectionForShell(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, G4int shellID, const G4Material *)
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:611
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4double bindingEnergy(G4int A, G4int Z)

◆ SelectStationary()

void G4DNARuddIonisationExtendedModel::SelectStationary ( G4bool  input)
inline

Definition at line 171 of file G4DNARuddIonisationExtendedModel.hh.

172{
173 statCode = input;
174}

Member Data Documentation

◆ fParticleChangeForGamma

G4ParticleChangeForGamma* G4DNARuddIonisationExtendedModel::fParticleChangeForGamma
protected

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