Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Attributes
G4DNARuddIonisationExtendedModel Class Reference

#include <G4DNARuddIonisationExtendedModel.hh>

Inheritance diagram for G4DNARuddIonisationExtendedModel:
G4VEmModel

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)
 
- 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 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 *, 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 *)
 
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)
 
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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
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
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
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 SetForceBuildTable (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

Protected Attributes

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

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 46 of file G4DNARuddIonisationExtendedModel.hh.

Constructor & Destructor Documentation

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

Definition at line 47 of file G4DNARuddIonisationExtendedModel.cc.

References python.hepunit::eV, fParticleChangeForGamma, G4cout, G4endl, python.hepunit::keV, python.hepunit::MeV, and G4VEmModel::SetDeexcitationFlag().

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

Definition at line 93 of file G4DNARuddIonisationExtendedModel.cc.

94 {
95  // Cross section
96 
97  std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
98  for (pos = tableData.begin(); pos != tableData.end(); ++pos)
99  {
100  G4DNACrossSectionDataSet* table = pos->second;
101  delete table;
102  }
103 
104  // The following removal is forbidden G4VEnergyLossModel takes care of deletion
105  // however coverity will signal this as an error
106  //if (fAtomDeexcitation) {delete fAtomDeexcitation;}
107 
108 }

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 374 of file G4DNARuddIonisationExtendedModel.cc.

References python.hepunit::cm, python.hepunit::eV, FatalException, G4DNACrossSectionDataSet::FindValue(), G4cout, G4endl, G4Exception(), G4Material::GetIndex(), G4DNAGenericIonsManager::GetIon(), G4ParticleDefinition::GetParticleName(), G4DNAGenericIonsManager::Instance(), and G4Proton::ProtonDefinition().

379 {
380  if (verboseLevel > 3)
381  G4cout << "Calling CrossSectionPerVolume() of G4DNARuddIonisationExtendedModel" << G4endl;
382 
383  // Calculate total cross section for model
384 
385  G4DNAGenericIonsManager *instance;
387 
388  if (
389  particleDefinition != G4Proton::ProtonDefinition()
390  &&
391  particleDefinition != instance->GetIon("hydrogen")
392  &&
393  particleDefinition != instance->GetIon("alpha++")
394  &&
395  particleDefinition != instance->GetIon("alpha+")
396  &&
397  particleDefinition != instance->GetIon("helium")
398  &&
399  particleDefinition != instance->GetIon("carbon")
400  &&
401  particleDefinition != instance->GetIon("nitrogen")
402  &&
403  particleDefinition != instance->GetIon("oxygen")
404  &&
405  particleDefinition != instance->GetIon("iron")
406  )
407 
408  return 0;
409 
410  G4double lowLim = 0;
411 
412  if ( particleDefinition == G4Proton::ProtonDefinition()
413  || particleDefinition == instance->GetIon("hydrogen")
414  )
415 
416  lowLim = lowEnergyLimitOfModelForA[1];
417 
418  else if ( particleDefinition == instance->GetIon("alpha++")
419  || particleDefinition == instance->GetIon("alpha+")
420  || particleDefinition == instance->GetIon("helium")
421  )
422 
423  lowLim = lowEnergyLimitOfModelForA[4];
424 
425  else lowLim = lowEnergyLimitOfModelForA[5];
426 
427  G4double highLim = 0;
428  G4double sigma=0;
429 
430 
431  G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
432 
433  if(waterDensity!= 0.0)
434 // if (material == nistwater || material->GetBaseMaterial() == nistwater)
435  {
436  const G4String& particleName = particleDefinition->GetParticleName();
437 
438  std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
439  pos2 = highEnergyLimit.find(particleName);
440 
441  if (pos2 != highEnergyLimit.end())
442  {
443  highLim = pos2->second;
444  }
445 
446  if (k <= highLim)
447  {
448 
449  //SI : XS must not be zero otherwise sampling of secondaries method ignored
450 
451  if (k < lowLim) k = lowLim;
452 
453  //
454 
455  std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
456  pos = tableData.find(particleName);
457 
458  if (pos != tableData.end())
459  {
460  G4DNACrossSectionDataSet* table = pos->second;
461  if (table != 0)
462  {
463  sigma = table->FindValue(k);
464  }
465  }
466  else
467  {
468  G4Exception("G4DNARuddIonisationExtendedModel::CrossSectionPerVolume","em0002",
469  FatalException,"Model not applicable to particle type.");
470  }
471 
472  } // if (k >= lowLim && k < highLim)
473 
474  if (verboseLevel > 2)
475  {
476  G4cout << "__________________________________" << G4endl;
477  G4cout << "°°° G4DNARuddIonisationExtendedModel - XS INFO START" << G4endl;
478  G4cout << "°°° Kinetic energy(eV)=" << k/eV << " particle : " << particleDefinition->GetParticleName() << G4endl;
479  G4cout << "°°° Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
480  G4cout << "°°° Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
481  // G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
482  G4cout << "°°° G4DNARuddIonisationExtendedModel - XS INFO END" << G4endl;
483 
484  }
485 
486  } // if (waterMaterial)
487 
488  return sigma*waterDensity;
489 // return sigma*material->GetAtomicNumDensityVector()[1];
490 
491 }
size_t GetIndex() const
Definition: G4Material.hh:260
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
G4GLOB_DLL std::ostream G4cout
static G4DNAGenericIonsManager * Instance(void)
virtual G4double FindValue(G4double e, G4int componentId=0) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4ParticleDefinition * GetIon(const G4String &name)
void G4DNARuddIonisationExtendedModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 112 of file G4DNARuddIonisationExtendedModel.cc.

References G4LossTableManager::AtomDeexcitation(), python.hepunit::eV, fParticleChangeForGamma, G4cout, G4endl, G4ParticleDefinition::GetAtomicMass(), G4DNAGenericIonsManager::GetIon(), G4Material::GetMaterial(), G4DNAMolecularMaterial::GetNumMolPerVolTableFor(), G4VEmModel::GetParticleChangeForGamma(), G4ParticleDefinition::GetParticleName(), G4VEmModel::HighEnergyLimit(), G4DNAGenericIonsManager::Instance(), G4DNAMolecularMaterial::Instance(), G4LossTableManager::Instance(), python.hepunit::keV, G4DNACrossSectionDataSet::LoadData(), G4VEmModel::LowEnergyLimit(), python.hepunit::m, python.hepunit::MeV, G4InuclParticleNames::proton, G4Proton::ProtonDefinition(), G4VEmModel::SetHighEnergyLimit(), and G4VEmModel::SetLowEnergyLimit().

114 {
115 
116  if (verboseLevel > 3)
117  G4cout << "Calling G4DNARuddIonisationExtendedModel::Initialise()" << G4endl;
118 
119  // Energy limits
120 
121  G4String fileProton("dna/sigma_ionisation_p_rudd");
122  G4String fileHydrogen("dna/sigma_ionisation_h_rudd");
123  G4String fileAlphaPlusPlus("dna/sigma_ionisation_alphaplusplus_rudd");
124  G4String fileAlphaPlus("dna/sigma_ionisation_alphaplus_rudd");
125  G4String fileHelium("dna/sigma_ionisation_he_rudd");
126  G4String fileCarbon("dna/sigma_ionisation_c_rudd");
127  G4String fileNitrogen("dna/sigma_ionisation_n_rudd");
128  G4String fileOxygen("dna/sigma_ionisation_o_rudd");
129  G4String fileIron("dna/sigma_ionisation_fe_rudd");
130 
131  G4DNAGenericIonsManager *instance;
134  G4ParticleDefinition* hydrogenDef = instance->GetIon("hydrogen");
135  G4ParticleDefinition* alphaPlusPlusDef = instance->GetIon("alpha++");
136  G4ParticleDefinition* alphaPlusDef = instance->GetIon("alpha+");
137  G4ParticleDefinition* heliumDef = instance->GetIon("helium");
138  G4ParticleDefinition* carbonDef = instance->GetIon("carbon");
139  G4ParticleDefinition* nitrogenDef = instance->GetIon("nitrogen");
140  G4ParticleDefinition* oxygenDef = instance->GetIon("oxygen");
141  G4ParticleDefinition* ironDef = instance->GetIon("iron");
142 
144  G4String hydrogen;
145  G4String alphaPlusPlus;
146  G4String alphaPlus;
147  G4String helium;
148  G4String carbon;
149  G4String nitrogen;
150  G4String oxygen;
151  G4String iron;
152 
153  G4double scaleFactor = 1 * m*m;
154 
155  // LIMITS AND DATA
156 
157  proton = protonDef->GetParticleName();
158  tableFile[proton] = fileProton;
159  lowEnergyLimit[proton] = lowEnergyLimitForA[1];
160  highEnergyLimit[proton] = 500. * keV;
161 
162  // Cross section
163 
165  eV,
166  scaleFactor );
167  tableProton->LoadData(fileProton);
168  tableData[proton] = tableProton;
169 
170  // **********************************************************************************************
171 
172  hydrogen = hydrogenDef->GetParticleName();
173  tableFile[hydrogen] = fileHydrogen;
174 
175  lowEnergyLimit[hydrogen] = lowEnergyLimitForA[1];
176  highEnergyLimit[hydrogen] = 100. * MeV;
177 
178  // Cross section
179 
181  eV,
182  scaleFactor );
183  tableHydrogen->LoadData(fileHydrogen);
184 
185  tableData[hydrogen] = tableHydrogen;
186 
187  // **********************************************************************************************
188 
189  alphaPlusPlus = alphaPlusPlusDef->GetParticleName();
190  tableFile[alphaPlusPlus] = fileAlphaPlusPlus;
191 
192  lowEnergyLimit[alphaPlusPlus] = lowEnergyLimitForA[4];
193  highEnergyLimit[alphaPlusPlus] = 400. * MeV;
194 
195  // Cross section
196 
198  eV,
199  scaleFactor );
200  tableAlphaPlusPlus->LoadData(fileAlphaPlusPlus);
201 
202  tableData[alphaPlusPlus] = tableAlphaPlusPlus;
203 
204  // **********************************************************************************************
205 
206  alphaPlus = alphaPlusDef->GetParticleName();
207  tableFile[alphaPlus] = fileAlphaPlus;
208 
209  lowEnergyLimit[alphaPlus] = lowEnergyLimitForA[4];
210  highEnergyLimit[alphaPlus] = 400. * MeV;
211 
212  // Cross section
213 
215  eV,
216  scaleFactor );
217  tableAlphaPlus->LoadData(fileAlphaPlus);
218  tableData[alphaPlus] = tableAlphaPlus;
219 
220  // **********************************************************************************************
221 
222  helium = heliumDef->GetParticleName();
223  tableFile[helium] = fileHelium;
224 
225  lowEnergyLimit[helium] = lowEnergyLimitForA[4];
226  highEnergyLimit[helium] = 400. * MeV;
227 
228  // Cross section
229 
231  eV,
232  scaleFactor );
233  tableHelium->LoadData(fileHelium);
234  tableData[helium] = tableHelium;
235 
236  // **********************************************************************************************
237 
238  carbon = carbonDef->GetParticleName();
239  tableFile[carbon] = fileCarbon;
240 
241  lowEnergyLimit[carbon] = lowEnergyLimitForA[5] * particle->GetAtomicMass();
242  highEnergyLimit[carbon] = 1e6* particle->GetAtomicMass() * MeV;
243 
244  // Cross section
245 
247  eV,
248  scaleFactor );
249  tableCarbon->LoadData(fileCarbon);
250  tableData[carbon] = tableCarbon;
251 
252  // **********************************************************************************************
253 
254  oxygen = oxygenDef->GetParticleName();
255  tableFile[oxygen] = fileOxygen;
256 
257  lowEnergyLimit[oxygen] = lowEnergyLimitForA[5]* particle->GetAtomicMass();
258  highEnergyLimit[oxygen] = 1e6* particle->GetAtomicMass()* MeV;
259 
260  // Cross section
261 
263  eV,
264  scaleFactor );
265  tableOxygen->LoadData(fileOxygen);
266  tableData[oxygen] = tableOxygen;
267 
268  // **********************************************************************************************
269 
270  nitrogen = nitrogenDef->GetParticleName();
271  tableFile[nitrogen] = fileNitrogen;
272 
273  lowEnergyLimit[nitrogen] = lowEnergyLimitForA[5]* particle->GetAtomicMass();
274  highEnergyLimit[nitrogen] = 1e6* particle->GetAtomicMass()* MeV;
275 
276  // Cross section
277 
279  eV,
280  scaleFactor );
281  tableNitrogen->LoadData(fileNitrogen);
282  tableData[nitrogen] = tableNitrogen;
283 
284  // **********************************************************************************************
285 
286  iron = ironDef->GetParticleName();
287  tableFile[iron] = fileIron;
288 
289  lowEnergyLimit[iron] = lowEnergyLimitForA[5]* particle->GetAtomicMass();
290  highEnergyLimit[iron] = 1e6* particle->GetAtomicMass()* MeV;
291 
292  // Cross section
293 
295  eV,
296  scaleFactor );
297  tableIron->LoadData(fileIron);
298  tableData[iron] = tableIron;
299 
300  // ZF Following lines can be replaced by:
301  SetLowEnergyLimit(lowEnergyLimit[particle->GetParticleName()]);
302  SetHighEnergyLimit(highEnergyLimit[particle->GetParticleName()]);
303  // at least for HZE
304 
305  /*
306  if (particle==protonDef)
307  {
308  SetLowEnergyLimit(lowEnergyLimit[proton]);
309  SetHighEnergyLimit(highEnergyLimit[proton]);
310  }
311 
312  if (particle==hydrogenDef)
313  {
314  SetLowEnergyLimit(lowEnergyLimit[hydrogen]);
315  SetHighEnergyLimit(highEnergyLimit[hydrogen]);
316  }
317 
318  if (particle==heliumDef)
319  {
320  SetLowEnergyLimit(lowEnergyLimit[helium]);
321  SetHighEnergyLimit(highEnergyLimit[helium]);
322  }
323 
324  if (particle==alphaPlusDef)
325  {
326  SetLowEnergyLimit(lowEnergyLimit[alphaPlus]);
327  SetHighEnergyLimit(highEnergyLimit[alphaPlus]);
328  }
329 
330  if (particle==alphaPlusPlusDef)
331  {
332  SetLowEnergyLimit(lowEnergyLimit[alphaPlusPlus]);
333  SetHighEnergyLimit(highEnergyLimit[alphaPlusPlus]);
334  }
335 
336  if (particle==carbonDef)
337  {
338  SetLowEnergyLimit(lowEnergyLimit[carbon]);
339  SetHighEnergyLimit(highEnergyLimit[carbon]);
340  }
341 
342  if (particle==oxygenDef)
343  {
344  SetLowEnergyLimit(lowEnergyLimit[oxygen]);
345  SetHighEnergyLimit(highEnergyLimit[oxygen]);
346  }*/
347 
348  //----------------------------------------------------------------------
349 
350  if( verboseLevel>0 )
351  {
352  G4cout << "Rudd ionisation model is initialized " << G4endl
353  << "Energy range: "
354  << LowEnergyLimit() / eV << " eV - "
355  << HighEnergyLimit() / keV << " keV for "
356  << particle->GetParticleName()
357  << G4endl;
358  }
359 
360  // Initialize water density pointer
362 
363  //
364 
365  fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
366 
367  if (isInitialised) { return; }
369  isInitialised = true;
370 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:599
static G4LossTableManager * Instance()
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:578
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
virtual G4bool LoadData(const G4String &argFileName)
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
G4GLOB_DLL std::ostream G4cout
const std::vector< double > * GetNumMolPerVolTableFor(const G4Material *) const
static G4DNAGenericIonsManager * Instance(void)
G4int GetAtomicMass() const
static G4DNAMolecularMaterial * Instance()
#define G4endl
Definition: G4ios.hh:61
G4VAtomDeexcitation * AtomDeexcitation()
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:690
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:121
G4ParticleDefinition * GetIon(const G4String &name)
void G4DNARuddIonisationExtendedModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple ,
const G4DynamicParticle particle,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 495 of file G4DNARuddIonisationExtendedModel.cc.

References G4InuclSpecialFunctions::bindingEnergy(), G4DNAChemistryManager::CreateWaterMolecule(), eIonizedMolecule, G4Electron::Electron(), fKShell, fParticleChangeForGamma, fStopAndKill, G4cout, G4endl, G4VAtomDeexcitation::GenerateParticles(), G4ParticleDefinition::GetAtomicMass(), G4VAtomDeexcitation::GetAtomicShell(), G4ParticleChangeForGamma::GetCurrentTrack(), G4DynamicParticle::GetDefinition(), G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), G4ParticleDefinition::GetParticleName(), G4DNAChemistryManager::Instance(), G4DNAWaterIonisationStructure::IonisationEnergy(), G4VParticleChange::ProposeLocalEnergyDeposit(), G4ParticleChangeForGamma::ProposeMomentumDirection(), G4VParticleChange::ProposeTrackStatus(), CLHEP::Hep3Vector::rotateUz(), and G4ParticleChangeForGamma::SetProposedKineticEnergy().

500 {
501  if (verboseLevel > 3)
502  G4cout << "Calling SampleSecondaries() of G4DNARuddIonisationExtendedModel" << G4endl;
503 
504  G4double lowLim = 0;
505  G4double highLim = 0;
506 
507  // ZF: the following line summarizes the commented part
508  if(particle->GetDefinition()->GetAtomicMass() <= 4) lowLim = killBelowEnergyForA[particle->GetDefinition()->GetAtomicMass()];
509  else lowLim = killBelowEnergyForA[5]*particle->GetDefinition()->GetAtomicMass();
510 
511  /*
512  if(particle->GetDefinition()->GetAtomicMass() >= 5) lowLim = killBelowEnergyForA[5]*particle->GetDefinition()->GetAtomicMass();
513 
514 
515  if ( particle->GetDefinition() == G4Proton::ProtonDefinition()
516  || particle->GetDefinition() == instance->GetIon("hydrogen")
517  )
518 
519  lowLim = killBelowEnergyForA[1];
520 
521  if ( particle->GetDefinition() == instance->GetIon("alpha++")
522  || particle->GetDefinition() == instance->GetIon("alpha+")
523  || particle->GetDefinition() == instance->GetIon("helium")
524  )
525 
526  lowLim = killBelowEnergyForA[4];*/
527 
528  G4double k = particle->GetKineticEnergy();
529 
530  const G4String& particleName = particle->GetDefinition()->GetParticleName();
531 
532  // SI - the following is useless since lowLim is already defined
533  /*
534  std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
535  pos1 = lowEnergyLimit.find(particleName);
536 
537  if (pos1 != lowEnergyLimit.end())
538  {
539  lowLim = pos1->second;
540  }
541  */
542 
543  std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
544  pos2 = highEnergyLimit.find(particleName);
545 
546  if (pos2 != highEnergyLimit.end())highLim = pos2->second;
547 
548  if (k >= lowLim && k < highLim)
549  {
550  G4ParticleDefinition* definition = particle->GetDefinition();
551  G4ParticleMomentum primaryDirection = particle->GetMomentumDirection();
552  /*
553  G4double particleMass = definition->GetPDGMass();
554  G4double totalEnergy = k + particleMass;
555  G4double pSquare = k*(totalEnergy+particleMass);
556  G4double totalMomentum = std::sqrt(pSquare);
557  */
558 
559  G4int ionizationShell = RandomSelect(k,particleName);
560 
561  // sample deexcitation
562  // here we assume that H_{2}O electronic levels are the same of Oxigen.
563  // this can be considered true with a rough 10% error in energy on K-shell,
564 
565  G4int secNumberInit = 0; // need to know at a certain point the enrgy of secondaries
566  G4int secNumberFinal = 0; // So I'll make the diference and then sum the energies
568  bindingEnergy = waterStructure.IonisationEnergy(ionizationShell);
569 
570  if(fAtomDeexcitation) {
571  G4int Z = 8;
573 
574  if (ionizationShell <5 && ionizationShell >1)
575  {
576  as = G4AtomicShellEnumerator(4-ionizationShell);
577  }
578  else if (ionizationShell <2)
579  {
580  as = G4AtomicShellEnumerator(3);
581  }
582 
583 
584  // DEBUG
585  // if (ionizationShell == 4) {
586  //
587  // G4cout << "Z: " << Z << " as: " << as
588  // << " ionizationShell: " << ionizationShell << " bindingEnergy: "<< bindingEnergy/eV << G4endl;
589  // G4cout << "Press <Enter> key to continue..." << G4endl;
590  // G4cin.ignore();
591  // }
592 
593 
594  const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
595  secNumberInit = fvect->size();
596  fAtomDeexcitation->GenerateParticles(fvect, shell, Z, 0, 0);
597  secNumberFinal = fvect->size();
598  }
599 
600 
601  G4double secondaryKinetic = RandomizeEjectedElectronEnergy(definition,k,ionizationShell);
602 
603  G4double cosTheta = 0.;
604  G4double phi = 0.;
605  RandomizeEjectedElectronDirection(definition, k,secondaryKinetic, cosTheta, phi, ionizationShell);
606 
607  G4double sinTheta = std::sqrt(1.-cosTheta*cosTheta);
608  G4double dirX = sinTheta*std::cos(phi);
609  G4double dirY = sinTheta*std::sin(phi);
610  G4double dirZ = cosTheta;
611  G4ThreeVector deltaDirection(dirX,dirY,dirZ);
612  deltaDirection.rotateUz(primaryDirection);
613 
614  // Ignored for ions on electrons
615  /*
616  G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 ));
617 
618  G4double finalPx = totalMomentum*primaryDirection.x() - deltaTotalMomentum*deltaDirection.x();
619  G4double finalPy = totalMomentum*primaryDirection.y() - deltaTotalMomentum*deltaDirection.y();
620  G4double finalPz = totalMomentum*primaryDirection.z() - deltaTotalMomentum*deltaDirection.z();
621  G4double finalMomentum = std::sqrt(finalPx*finalPx+finalPy*finalPy+finalPz*finalPz);
622  finalPx /= finalMomentum;
623  finalPy /= finalMomentum;
624  finalPz /= finalMomentum;
625 
626  G4ThreeVector direction;
627  direction.set(finalPx,finalPy,finalPz);
628 
629  fParticleChangeForGamma->ProposeMomentumDirection(direction.unit()) ;
630  */
631 
633  G4double scatteredEnergy = k-bindingEnergy-secondaryKinetic;
634  G4double deexSecEnergy = 0;
635  for (G4int j=secNumberInit; j < secNumberFinal; j++) {
636 
637  deexSecEnergy = deexSecEnergy + (*fvect)[j]->GetKineticEnergy();
638 
639  }
640 
642  fParticleChangeForGamma->ProposeLocalEnergyDeposit(k-scatteredEnergy-secondaryKinetic-deexSecEnergy);
643 
644  G4DynamicParticle* dp = new G4DynamicParticle (G4Electron::Electron(),deltaDirection,secondaryKinetic) ;
645  fvect->push_back(dp);
646 
647  const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
649  ionizationShell,
650  theIncomingTrack);
651  }
652 
653  // SI - not useful since low energy of model is 0 eV
654 
655  if (k < lowLim)
656  {
660  }
661 
662 }
G4double GetKineticEnergy() const
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
const G4String & GetParticleName() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
G4int GetAtomicMass() const
static G4DNAChemistryManager * Instance()
void CreateWaterMolecule(ElectronicModification, G4int, const G4Track *)
const G4Track * GetCurrentTrack() const
static G4Electron * Electron()
Definition: G4Electron.cc:94
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
G4double bindingEnergy(G4int A, G4int Z)
G4AtomicShellEnumerator

Field Documentation

G4ParticleChangeForGamma* G4DNARuddIonisationExtendedModel::fParticleChangeForGamma
protected

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