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

#include <G4DNARuddIonisationModel.hh>

Inheritance diagram for G4DNARuddIonisationModel:
G4VEmModel

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)
 
- 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 G4DNARuddIonisationModel.hh.

Constructor & Destructor Documentation

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

Definition at line 44 of file G4DNARuddIonisationModel.cc.

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

46  :G4VEmModel(nam),isInitialised(false)
47 {
48  // nistwater = G4NistManager::Instance()->FindOrBuildMaterial("G4_WATER");
49  fpWaterDensity = 0;
50 
51  slaterEffectiveCharge[0]=0.;
52  slaterEffectiveCharge[1]=0.;
53  slaterEffectiveCharge[2]=0.;
54  sCoefficient[0]=0.;
55  sCoefficient[1]=0.;
56  sCoefficient[2]=0.;
57 
58  lowEnergyLimitForZ1 = 0 * eV;
59  lowEnergyLimitForZ2 = 0 * eV;
60  lowEnergyLimitOfModelForZ1 = 100 * eV;
61  lowEnergyLimitOfModelForZ2 = 1 * keV;
62  killBelowEnergyForZ1 = lowEnergyLimitOfModelForZ1;
63  killBelowEnergyForZ2 = lowEnergyLimitOfModelForZ2;
64 
65  verboseLevel= 0;
66  // Verbosity scale:
67  // 0 = nothing
68  // 1 = warning for energy non-conservation
69  // 2 = details of energy budget
70  // 3 = calculation of cross sections, file openings, sampling of atoms
71  // 4 = entering in methods
72 
73  if( verboseLevel>0 )
74  {
75  G4cout << "Rudd ionisation model is constructed " << G4endl;
76  }
77 
78  //Mark this model as "applicable" for atomic deexcitation
79  SetDeexcitationFlag(true);
80  fAtomDeexcitation = 0;
82 }
G4ParticleChangeForGamma * fParticleChangeForGamma
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
G4DNARuddIonisationModel::~G4DNARuddIonisationModel ( )
virtual

Definition at line 86 of file G4DNARuddIonisationModel.cc.

87 {
88  // Cross section
89 
90  std::map< G4String,G4DNACrossSectionDataSet*,std::less<G4String> >::iterator pos;
91  for (pos = tableData.begin(); pos != tableData.end(); ++pos)
92  {
93  G4DNACrossSectionDataSet* table = pos->second;
94  delete table;
95  }
96 
97 
98  // The following removal is forbidden since G4VEnergyLossmodel takes care of deletion
99  // Coverity however will signal this as an error
100  //if (fAtomDeexcitation) {delete fAtomDeexcitation;}
101 
102 }

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 278 of file G4DNARuddIonisationModel.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().

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

Implements G4VEmModel.

Definition at line 106 of file G4DNARuddIonisationModel.cc.

References G4LossTableManager::AtomDeexcitation(), python.hepunit::eV, fParticleChangeForGamma, G4cout, G4endl, 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().

108 {
109 
110  if (verboseLevel > 3)
111  G4cout << "Calling G4DNARuddIonisationModel::Initialise()" << G4endl;
112 
113  // Energy limits
114 
115  G4String fileProton("dna/sigma_ionisation_p_rudd");
116  G4String fileHydrogen("dna/sigma_ionisation_h_rudd");
117  G4String fileAlphaPlusPlus("dna/sigma_ionisation_alphaplusplus_rudd");
118  G4String fileAlphaPlus("dna/sigma_ionisation_alphaplus_rudd");
119  G4String fileHelium("dna/sigma_ionisation_he_rudd");
120 
121  G4DNAGenericIonsManager *instance;
124  G4ParticleDefinition* hydrogenDef = instance->GetIon("hydrogen");
125  G4ParticleDefinition* alphaPlusPlusDef = instance->GetIon("alpha++");
126  G4ParticleDefinition* alphaPlusDef = instance->GetIon("alpha+");
127  G4ParticleDefinition* heliumDef = instance->GetIon("helium");
128 
130  G4String hydrogen;
131  G4String alphaPlusPlus;
132  G4String alphaPlus;
133  G4String helium;
134 
135  G4double scaleFactor = 1 * m*m;
136 
137  // LIMITS AND DATA
138 
139  // ********************************************************
140 
141  proton = protonDef->GetParticleName();
142  tableFile[proton] = fileProton;
143 
144  lowEnergyLimit[proton] = lowEnergyLimitForZ1;
145  highEnergyLimit[proton] = 500. * keV;
146 
147  // Cross section
148 
150  eV,
151  scaleFactor );
152  tableProton->LoadData(fileProton);
153  tableData[proton] = tableProton;
154 
155  // ********************************************************
156 
157  hydrogen = hydrogenDef->GetParticleName();
158  tableFile[hydrogen] = fileHydrogen;
159 
160  lowEnergyLimit[hydrogen] = lowEnergyLimitForZ1;
161  highEnergyLimit[hydrogen] = 100. * MeV;
162 
163  // Cross section
164 
166  eV,
167  scaleFactor );
168  tableHydrogen->LoadData(fileHydrogen);
169 
170  tableData[hydrogen] = tableHydrogen;
171 
172  // ********************************************************
173 
174  alphaPlusPlus = alphaPlusPlusDef->GetParticleName();
175  tableFile[alphaPlusPlus] = fileAlphaPlusPlus;
176 
177  lowEnergyLimit[alphaPlusPlus] = lowEnergyLimitForZ2;
178  highEnergyLimit[alphaPlusPlus] = 400. * MeV;
179 
180  // Cross section
181 
183  eV,
184  scaleFactor );
185  tableAlphaPlusPlus->LoadData(fileAlphaPlusPlus);
186 
187  tableData[alphaPlusPlus] = tableAlphaPlusPlus;
188 
189  // ********************************************************
190 
191  alphaPlus = alphaPlusDef->GetParticleName();
192  tableFile[alphaPlus] = fileAlphaPlus;
193 
194  lowEnergyLimit[alphaPlus] = lowEnergyLimitForZ2;
195  highEnergyLimit[alphaPlus] = 400. * MeV;
196 
197  // Cross section
198 
200  eV,
201  scaleFactor );
202  tableAlphaPlus->LoadData(fileAlphaPlus);
203  tableData[alphaPlus] = tableAlphaPlus;
204 
205  // ********************************************************
206 
207  helium = heliumDef->GetParticleName();
208  tableFile[helium] = fileHelium;
209 
210  lowEnergyLimit[helium] = lowEnergyLimitForZ2;
211  highEnergyLimit[helium] = 400. * MeV;
212 
213  // Cross section
214 
216  eV,
217  scaleFactor );
218  tableHelium->LoadData(fileHelium);
219  tableData[helium] = tableHelium;
220 
221  //
222 
223  if (particle==protonDef)
224  {
225  SetLowEnergyLimit(lowEnergyLimit[proton]);
226  SetHighEnergyLimit(highEnergyLimit[proton]);
227  }
228 
229  if (particle==hydrogenDef)
230  {
231  SetLowEnergyLimit(lowEnergyLimit[hydrogen]);
232  SetHighEnergyLimit(highEnergyLimit[hydrogen]);
233  }
234 
235  if (particle==heliumDef)
236  {
237  SetLowEnergyLimit(lowEnergyLimit[helium]);
238  SetHighEnergyLimit(highEnergyLimit[helium]);
239  }
240 
241  if (particle==alphaPlusDef)
242  {
243  SetLowEnergyLimit(lowEnergyLimit[alphaPlus]);
244  SetHighEnergyLimit(highEnergyLimit[alphaPlus]);
245  }
246 
247  if (particle==alphaPlusPlusDef)
248  {
249  SetLowEnergyLimit(lowEnergyLimit[alphaPlusPlus]);
250  SetHighEnergyLimit(highEnergyLimit[alphaPlusPlus]);
251  }
252 
253  if( verboseLevel>0 )
254  {
255  G4cout << "Rudd ionisation model is initialized " << G4endl
256  << "Energy range: "
257  << LowEnergyLimit() / eV << " eV - "
258  << HighEnergyLimit() / keV << " keV for "
259  << particle->GetParticleName()
260  << G4endl;
261  }
262 
263  // Initialize water density pointer
265 
266  //
267 
268  fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
269 
270  if (isInitialised) { return; }
272  isInitialised = true;
273 
274 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:599
static G4LossTableManager * Instance()
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
G4ParticleChangeForGamma * fParticleChangeForGamma
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)
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 G4DNARuddIonisationModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple ,
const G4DynamicParticle particle,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 404 of file G4DNARuddIonisationModel.cc.

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

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

Field Documentation

G4ParticleChangeForGamma* G4DNARuddIonisationModel::fParticleChangeForGamma
protected

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