Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4DNARuddIonisationModel Class Reference

#include <G4DNARuddIonisationModel.hh>

Inheritance diagram for G4DNARuddIonisationModel:
Collaboration diagram for G4DNARuddIonisationModel:

Public Member Functions

 G4DNARuddIonisationModel (const G4ParticleDefinition *p=0, const G4String &nam="DNARuddIonisationModel")
 
virtual ~G4DNARuddIonisationModel ()
 
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 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 *, 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)
 
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=nullptr)
 
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
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 

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 *)
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 46 of file G4DNARuddIonisationModel.hh.

Constructor & Destructor Documentation

G4DNARuddIonisationModel::G4DNARuddIonisationModel ( const G4ParticleDefinition p = 0,
const G4String nam = "DNARuddIonisationModel" 
)

Definition at line 47 of file G4DNARuddIonisationModel.cc.

48  :
49  G4VEmModel(nam), isInitialised(false)
50 {
51  fpWaterDensity = 0;
52 
53  slaterEffectiveCharge[0] = 0.;
54  slaterEffectiveCharge[1] = 0.;
55  slaterEffectiveCharge[2] = 0.;
56  sCoefficient[0] = 0.;
57  sCoefficient[1] = 0.;
58  sCoefficient[2] = 0.;
59 
60  lowEnergyLimitForZ1 = 0 * eV;
61  lowEnergyLimitForZ2 = 0 * eV;
62  lowEnergyLimitOfModelForZ1 = 100 * eV;
63  lowEnergyLimitOfModelForZ2 = 1 * keV;
64  killBelowEnergyForZ1 = lowEnergyLimitOfModelForZ1;
65  killBelowEnergyForZ2 = lowEnergyLimitOfModelForZ2;
66 
67  verboseLevel = 0;
68  // Verbosity scale:
69  // 0 = nothing
70  // 1 = warning for energy non-conservation
71  // 2 = details of energy budget
72  // 3 = calculation of cross sections, file openings, sampling of atoms
73  // 4 = entering in methods
74 
75  if (verboseLevel > 0)
76  {
77  G4cout << "Rudd ionisation model is constructed " << G4endl;
78  }
79 
80  // define default angular generator
82 
83  //Mark this model as "applicable" for atomic deexcitation
84  SetDeexcitationFlag(true);
85  fAtomDeexcitation = 0;
87 
88  // Selection of stationary mode
89 
90  statCode = false;
91 }
G4ParticleChangeForGamma * fParticleChangeForGamma
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
G4GLOB_DLL std::ostream G4cout
static constexpr double eV
Definition: G4SIunits.hh:215
void SetAngularDistribution(G4VEmAngularDistribution *)
Definition: G4VEmModel.hh:623
#define G4endl
Definition: G4ios.hh:61
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:780
static constexpr double keV
Definition: G4SIunits.hh:216

Here is the call graph for this function:

G4DNARuddIonisationModel::~G4DNARuddIonisationModel ( )
virtual

Definition at line 95 of file G4DNARuddIonisationModel.cc.

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

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 289 of file G4DNARuddIonisationModel.cc.

294 {
295  if (verboseLevel > 3)
296  {
297  G4cout << "Calling CrossSectionPerVolume() of G4DNARuddIonisationModel"
298  << G4endl;
299  }
300 
301  // Calculate total cross section for model
302 
303  G4DNAGenericIonsManager *instance;
305 
306  if (
307  particleDefinition != G4Proton::ProtonDefinition()
308  &&
309  particleDefinition != instance->GetIon("hydrogen")
310  &&
311  particleDefinition != instance->GetIon("alpha++")
312  &&
313  particleDefinition != instance->GetIon("alpha+")
314  &&
315  particleDefinition != instance->GetIon("helium")
316  )
317 
318  return 0;
319 
320  G4double lowLim = 0;
321 
322  if ( particleDefinition == G4Proton::ProtonDefinition()
323  || particleDefinition == instance->GetIon("hydrogen")
324  )
325 
326  lowLim = lowEnergyLimitOfModelForZ1;
327 
328  if ( particleDefinition == instance->GetIon("alpha++")
329  || particleDefinition == instance->GetIon("alpha+")
330  || particleDefinition == instance->GetIon("helium")
331  )
332 
333  lowLim = lowEnergyLimitOfModelForZ2;
334 
335  G4double highLim = 0;
336  G4double sigma=0;
337 
338  G4double waterDensity = (*fpWaterDensity)[material->GetIndex()];
339 
340  if(waterDensity!= 0.0)
341  // if (material == nistwater || material->GetBaseMaterial() == nistwater)
342  {
343  const G4String& particleName = particleDefinition->GetParticleName();
344 
345  // SI - the following is useless since lowLim is already defined
346  /*
347  std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
348  pos1 = lowEnergyLimit.find(particleName);
349 
350  if (pos1 != lowEnergyLimit.end())
351  {
352  lowLim = pos1->second;
353  }
354  */
355 
356  std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
357  pos2 = highEnergyLimit.find(particleName);
358 
359  if (pos2 != highEnergyLimit.end())
360  {
361  highLim = pos2->second;
362  }
363 
364  if (k <= highLim)
365  {
366  //SI : XS must not be zero otherwise sampling of secondaries method ignored
367 
368  if (k < lowLim) k = lowLim;
369 
370  //
371 
372  std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
373  pos = tableData.find(particleName);
374 
375  if (pos != tableData.end())
376  {
377  G4DNACrossSectionDataSet* table = pos->second;
378  if (table != 0)
379  {
380  sigma = table->FindValue(k);
381  }
382  }
383  else
384  {
385  G4Exception("G4DNARuddIonisationModel::CrossSectionPerVolume","em0002",
386  FatalException,"Model not applicable to particle type.");
387  }
388 
389  } // if (k >= lowLim && k < highLim)
390 
391  if (verboseLevel > 2)
392  {
393  G4cout << "__________________________________" << G4endl;
394  G4cout << "G4DNARuddIonisationModel - XS INFO START" << G4endl;
395  G4cout << "Kinetic energy(eV)=" << k/eV << " particle : " << particleDefinition->GetParticleName() << G4endl;
396  G4cout << "Cross section per water molecule (cm^2)=" << sigma/cm/cm << G4endl;
397  G4cout << "Cross section per water molecule (cm^-1)=" << sigma*waterDensity/(1./cm) << G4endl;
398  // G4cout << " - Cross section per water molecule (cm^-1)=" << sigma*material->GetAtomicNumDensityVector()[1]/(1./cm) << G4endl;
399  G4cout << "G4DNARuddIonisationModel - XS INFO END" << G4endl;
400  }
401 
402  } // if (waterMaterial)
403  else
404  {
405  if (verboseLevel > 2)
406  {
407  G4cout << "Warning : RuddIonisationModel: WATER DENSITY IS NULL" << G4endl;
408  }
409  }
410 
411  return sigma*waterDensity;
412  // return sigma*material->GetAtomicNumDensityVector()[1];
413 
414 }
size_t GetIndex() const
Definition: G4Material.hh:262
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
G4GLOB_DLL std::ostream G4cout
static constexpr double cm
Definition: G4SIunits.hh:119
static constexpr double eV
Definition: G4SIunits.hh:215
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
static const G4double pos
G4ParticleDefinition * GetIon(const G4String &name)

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 114 of file G4DNARuddIonisationModel.cc.

116 {
117 
118  if (verboseLevel > 3)
119  {
120  G4cout << "Calling G4DNARuddIonisationModel::Initialise()" << G4endl;
121  }
122 
123  // Energy limits
124 
125  G4String fileProton("dna/sigma_ionisation_p_rudd");
126  G4String fileHydrogen("dna/sigma_ionisation_h_rudd");
127  G4String fileAlphaPlusPlus("dna/sigma_ionisation_alphaplusplus_rudd");
128  G4String fileAlphaPlus("dna/sigma_ionisation_alphaplus_rudd");
129  G4String fileHelium("dna/sigma_ionisation_he_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 
140  G4String hydrogen;
141  G4String alphaPlusPlus;
142  G4String alphaPlus;
143  G4String helium;
144 
145  G4double scaleFactor = 1 * m*m;
146 
147  // LIMITS AND DATA
148 
149  // ********************************************************
150 
151  proton = protonDef->GetParticleName();
152  tableFile[proton] = fileProton;
153 
154  lowEnergyLimit[proton] = lowEnergyLimitForZ1;
155  highEnergyLimit[proton] = 500. * keV;
156 
157  // Cross section
158 
160  eV,
161  scaleFactor );
162  tableProton->LoadData(fileProton);
163  tableData[proton] = tableProton;
164 
165  // ********************************************************
166 
167  hydrogen = hydrogenDef->GetParticleName();
168  tableFile[hydrogen] = fileHydrogen;
169 
170  lowEnergyLimit[hydrogen] = lowEnergyLimitForZ1;
171  highEnergyLimit[hydrogen] = 100. * MeV;
172 
173  // Cross section
174 
176  eV,
177  scaleFactor );
178  tableHydrogen->LoadData(fileHydrogen);
179 
180  tableData[hydrogen] = tableHydrogen;
181 
182  // ********************************************************
183 
184  alphaPlusPlus = alphaPlusPlusDef->GetParticleName();
185  tableFile[alphaPlusPlus] = fileAlphaPlusPlus;
186 
187  lowEnergyLimit[alphaPlusPlus] = lowEnergyLimitForZ2;
188  highEnergyLimit[alphaPlusPlus] = 400. * MeV;
189 
190  // Cross section
191 
193  eV,
194  scaleFactor );
195  tableAlphaPlusPlus->LoadData(fileAlphaPlusPlus);
196 
197  tableData[alphaPlusPlus] = tableAlphaPlusPlus;
198 
199  // ********************************************************
200 
201  alphaPlus = alphaPlusDef->GetParticleName();
202  tableFile[alphaPlus] = fileAlphaPlus;
203 
204  lowEnergyLimit[alphaPlus] = lowEnergyLimitForZ2;
205  highEnergyLimit[alphaPlus] = 400. * MeV;
206 
207  // Cross section
208 
210  eV,
211  scaleFactor );
212  tableAlphaPlus->LoadData(fileAlphaPlus);
213  tableData[alphaPlus] = tableAlphaPlus;
214 
215  // ********************************************************
216 
217  helium = heliumDef->GetParticleName();
218  tableFile[helium] = fileHelium;
219 
220  lowEnergyLimit[helium] = lowEnergyLimitForZ2;
221  highEnergyLimit[helium] = 400. * MeV;
222 
223  // Cross section
224 
226  eV,
227  scaleFactor );
228  tableHelium->LoadData(fileHelium);
229  tableData[helium] = tableHelium;
230 
231  //
232 
233  if (particle==protonDef)
234  {
235  SetLowEnergyLimit(lowEnergyLimit[proton]);
236  SetHighEnergyLimit(highEnergyLimit[proton]);
237  }
238 
239  if (particle==hydrogenDef)
240  {
241  SetLowEnergyLimit(lowEnergyLimit[hydrogen]);
242  SetHighEnergyLimit(highEnergyLimit[hydrogen]);
243  }
244 
245  if (particle==heliumDef)
246  {
247  SetLowEnergyLimit(lowEnergyLimit[helium]);
248  SetHighEnergyLimit(highEnergyLimit[helium]);
249  }
250 
251  if (particle==alphaPlusDef)
252  {
253  SetLowEnergyLimit(lowEnergyLimit[alphaPlus]);
254  SetHighEnergyLimit(highEnergyLimit[alphaPlus]);
255  }
256 
257  if (particle==alphaPlusPlusDef)
258  {
259  SetLowEnergyLimit(lowEnergyLimit[alphaPlusPlus]);
260  SetHighEnergyLimit(highEnergyLimit[alphaPlusPlus]);
261  }
262 
263  if( verboseLevel>0 )
264  {
265  G4cout << "Rudd ionisation model is initialized " << G4endl
266  << "Energy range: "
267  << LowEnergyLimit() / eV << " eV - "
268  << HighEnergyLimit() / keV << " keV for "
269  << particle->GetParticleName()
270  << G4endl;
271  }
272 
273  // Initialize water density pointer
275 
276  //
277 
278  fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
279 
280  if (isInitialised)
281  { return;}
283  isInitialised = true;
284 
285 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:640
static G4LossTableManager * Instance()
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:633
G4ParticleChangeForGamma * fParticleChangeForGamma
static G4Material * GetMaterial(const G4String &name, G4bool warning=true)
Definition: G4Material.cc:602
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
virtual G4bool LoadData(const G4String &argFileName)
const G4String & GetParticleName() const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:724
G4GLOB_DLL std::ostream G4cout
const std::vector< double > * GetNumMolPerVolTableFor(const G4Material *) const
static constexpr double m
Definition: G4SIunits.hh:129
static constexpr double eV
Definition: G4SIunits.hh:215
static G4DNAGenericIonsManager * Instance(void)
static G4DNAMolecularMaterial * Instance()
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
G4VAtomDeexcitation * AtomDeexcitation()
double G4double
Definition: G4Types.hh:76
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:731
static constexpr double keV
Definition: G4SIunits.hh:216
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132
G4ParticleDefinition * GetIon(const G4String &name)

Here is the call graph for this function:

void G4DNARuddIonisationModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle particle,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 418 of file G4DNARuddIonisationModel.cc.

423 {
424  if (verboseLevel > 3)
425  {
426  G4cout << "Calling SampleSecondaries() of G4DNARuddIonisationModel"
427  << G4endl;
428  }
429 
430  G4double lowLim = 0;
431  G4double highLim = 0;
432 
433  G4DNAGenericIonsManager *instance;
435 
436  if ( particle->GetDefinition() == G4Proton::ProtonDefinition()
437  || particle->GetDefinition() == instance->GetIon("hydrogen")
438  )
439 
440  lowLim = killBelowEnergyForZ1;
441 
442  if ( particle->GetDefinition() == instance->GetIon("alpha++")
443  || particle->GetDefinition() == instance->GetIon("alpha+")
444  || particle->GetDefinition() == instance->GetIon("helium")
445  )
446 
447  lowLim = killBelowEnergyForZ2;
448 
449  G4double k = particle->GetKineticEnergy();
450 
451  const G4String& particleName = particle->GetDefinition()->GetParticleName();
452 
453  // SI - the following is useless since lowLim is already defined
454  /*
455  std::map< G4String,G4double,std::less<G4String> >::iterator pos1;
456  pos1 = lowEnergyLimit.find(particleName);
457 
458  if (pos1 != lowEnergyLimit.end())
459  {
460  lowLim = pos1->second;
461  }
462  */
463 
464  std::map< G4String,G4double,std::less<G4String> >::iterator pos2;
465  pos2 = highEnergyLimit.find(particleName);
466 
467  if (pos2 != highEnergyLimit.end())
468  {
469  highLim = pos2->second;
470  }
471 
472  if (k >= lowLim && k <= highLim)
473  {
474  G4ParticleDefinition* definition = particle->GetDefinition();
475  G4ParticleMomentum primaryDirection = particle->GetMomentumDirection();
476  /*
477  G4double particleMass = definition->GetPDGMass();
478  G4double totalEnergy = k + particleMass;
479  G4double pSquare = k*(totalEnergy+particleMass);
480  G4double totalMomentum = std::sqrt(pSquare);
481  */
482 
483  G4int ionizationShell = RandomSelect(k,particleName);
484 
485  // sample deexcitation
486  // here we assume that H_{2}O electronic levels are the same of Oxigen.
487  // this can be considered true with a rough 10% error in energy on K-shell,
488 
489  G4int secNumberInit = 0;// need to know at a certain point the enrgy of secondaries
490  G4int secNumberFinal = 0;// So I'll make the diference and then sum the energies
491 
493  bindingEnergy = waterStructure.IonisationEnergy(ionizationShell);
494 
495  //SI: additional protection if tcs interpolation method is modified
496  if (k<bindingEnergy) return;
497  //
498 
499  G4int Z = 8;
500  if(fAtomDeexcitation)
501  {
503 
504  if (ionizationShell <5 && ionizationShell >1)
505  {
506  as = G4AtomicShellEnumerator(4-ionizationShell);
507  }
508  else if (ionizationShell <2)
509  {
510  as = G4AtomicShellEnumerator(3);
511  }
512 
513  // DEBUG
514  // if (ionizationShell == 4) {
515  //
516  // G4cout << "Z: " << Z << " as: " << as
517  // << " ionizationShell: " << ionizationShell << " bindingEnergy: "<< bindingEnergy/eV << G4endl;
518  // G4cout << "Press <Enter> key to continue..." << G4endl;
519  // G4cin.ignore();
520  // }
521 
522  const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
523  secNumberInit = fvect->size();
524  fAtomDeexcitation->GenerateParticles(fvect, shell, Z, 0, 0);
525  secNumberFinal = fvect->size();
526  }
527 
528  G4double secondaryKinetic = RandomizeEjectedElectronEnergy(definition,k,ionizationShell);
529 
530  G4ThreeVector deltaDirection =
531  GetAngularDistribution()->SampleDirectionForShell(particle, secondaryKinetic,
532  Z, ionizationShell,
533  couple->GetMaterial());
534 
535  // Ignored for ions on electrons
536  /*
537  G4double deltaTotalMomentum = std::sqrt(secondaryKinetic*(secondaryKinetic + 2.*electron_mass_c2 ));
538 
539  G4double finalPx = totalMomentum*primaryDirection.x() - deltaTotalMomentum*deltaDirection.x();
540  G4double finalPy = totalMomentum*primaryDirection.y() - deltaTotalMomentum*deltaDirection.y();
541  G4double finalPz = totalMomentum*primaryDirection.z() - deltaTotalMomentum*deltaDirection.z();
542  G4double finalMomentum = std::sqrt(finalPx*finalPx+finalPy*finalPy+finalPz*finalPz);
543  finalPx /= finalMomentum;
544  finalPy /= finalMomentum;
545  finalPz /= finalMomentum;
546 
547  G4ThreeVector direction;
548  direction.set(finalPx,finalPy,finalPz);
549 
550  fParticleChangeForGamma->ProposeMomentumDirection(direction.unit()) ;
551  */
553 
554  G4double scatteredEnergy = k-bindingEnergy-secondaryKinetic;
555  G4double deexSecEnergy = 0;
556  for (G4int j=secNumberInit; j < secNumberFinal; j++)
557  {
558 
559  deexSecEnergy = deexSecEnergy + (*fvect)[j]->GetKineticEnergy();
560 
561  }
562 
563  if (!statCode)
564  {
566  fParticleChangeForGamma->ProposeLocalEnergyDeposit(k-scatteredEnergy-secondaryKinetic-deexSecEnergy);
567  }
568  else
569  {
572  }
573 
574  // debug
575  // k-scatteredEnergy-secondaryKinetic-deexSecEnergy = k-(k-bindingEnergy-secondaryKinetic)-secondaryKinetic-deexSecEnergy =
576  // = k-k+bindingEnergy+secondaryKinetic-secondaryKinetic-deexSecEnergy=
577  // = bindingEnergy-deexSecEnergy
578  // SO deexSecEnergy=0 => LocalEnergyDeposit = bindingEnergy
579 
580  G4DynamicParticle* dp = new G4DynamicParticle (G4Electron::Electron(),deltaDirection,secondaryKinetic);
581  fvect->push_back(dp);
582 
583  const G4Track * theIncomingTrack = fParticleChangeForGamma->GetCurrentTrack();
585  ionizationShell,
586  theIncomingTrack);
587  }
588 
589  // SI - not useful since low energy of model is 0 eV
590 
591  if (k < lowLim)
592  {
596  }
597 
598 }
G4double GetKineticEnergy() const
G4ParticleChangeForGamma * fParticleChangeForGamma
static G4Proton * ProtonDefinition()
Definition: G4Proton.cc:88
G4VEmAngularDistribution * GetAngularDistribution()
Definition: G4VEmModel.hh:616
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 G4ThreeVector & SampleDirectionForShell(const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, G4int shellID, const G4Material *)
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
static G4DNAGenericIonsManager * Instance(void)
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
const G4Material * GetMaterial() const
G4ParticleDefinition * GetIon(const G4String &name)

Here is the call graph for this function:

void G4DNARuddIonisationModel::SelectStationary ( G4bool  input)
inline

Definition at line 163 of file G4DNARuddIonisationModel.hh.

164 {
165  statCode = input;
166 }

Member Data Documentation

G4ParticleChangeForGamma* G4DNARuddIonisationModel::fParticleChangeForGamma
protected

Definition at line 74 of file G4DNARuddIonisationModel.hh.


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