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

#include <G4LivermoreIonisationModel.hh>

Inheritance diagram for G4LivermoreIonisationModel:
G4VEmModel

Public Member Functions

 G4LivermoreIonisationModel (const G4ParticleDefinition *p=0, const G4String &processName="LowEnergyIoni")
 
virtual ~G4LivermoreIonisationModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
 
void SetVerboseLevel (G4int vl)
 
G4int GetVerboseLevel ()
 
- 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 CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, 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

G4ParticleChangeForLossfParticleChange
 
- 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 58 of file G4LivermoreIonisationModel.hh.

Constructor & Destructor Documentation

G4LivermoreIonisationModel::G4LivermoreIonisationModel ( const G4ParticleDefinition p = 0,
const G4String processName = "LowEnergyIoni" 
)

Definition at line 77 of file G4LivermoreIonisationModel.cc.

References python.hepunit::eV, python.hepunit::GeV, and G4AtomicTransitionManager::Instance().

78  :
79  G4VEmModel(nam), fParticleChange(0),
80  isInitialised(false),
81  crossSectionHandler(0), energySpectrum(0)
82 {
83  fIntrinsicLowEnergyLimit = 10.0*eV;
84  fIntrinsicHighEnergyLimit = 100.0*GeV;
85 
86  verboseLevel = 0;
87  //SetAngularDistribution(new G4DeltaAngle());
88 
89  transitionManager = G4AtomicTransitionManager::Instance();
90 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:65
G4ParticleChangeForLoss * fParticleChange
static G4AtomicTransitionManager * Instance()
G4LivermoreIonisationModel::~G4LivermoreIonisationModel ( )
virtual

Definition at line 94 of file G4LivermoreIonisationModel.cc.

95 {
96  delete energySpectrum;
97  delete crossSectionHandler;
98 }

Member Function Documentation

G4double G4LivermoreIonisationModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  kinEnergy,
G4double  Z,
G4double  A = 0,
G4double  cut = 0,
G4double  emax = DBL_MAX 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 174 of file G4LivermoreIonisationModel.cc.

References python.hepunit::barn, FatalException, G4cout, G4endl, G4Exception(), G4eIonisationCrossSectionHandler::GetCrossSectionAboveThresholdForElement(), and python.hepunit::keV.

180 {
181  G4int iZ = (G4int) Z;
182  if (!crossSectionHandler)
183  {
184  G4Exception("G4LivermoreIonisationModel::ComputeCrossSectionPerAtom",
185  "em1007",FatalException,
186  "The cross section handler is not correctly initialized");
187  return 0;
188  }
189 
190  //The cut is already included in the crossSectionHandler
191  G4double cs =
192  crossSectionHandler->GetCrossSectionAboveThresholdForElement(energy,
193  cutEnergy,
194  iZ);
195 
196  if (verboseLevel > 1)
197  {
198  G4cout << "G4LivermoreIonisationModel " << G4endl;
199  G4cout << "Cross section for delta emission > "
200  << cutEnergy/keV << " keV at "
201  << energy/keV << " keV and Z = " << iZ << " --> "
202  << cs/barn << " barn" << G4endl;
203  }
204  return cs;
205 }
G4double GetCrossSectionAboveThresholdForElement(G4double energy, G4double cutEnergy, G4int Z)
int G4int
Definition: G4Types.hh:78
double precision function energy(A, Z)
Definition: dpm25nuc6.f:4106
G4GLOB_DLL std::ostream G4cout
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
G4double G4LivermoreIonisationModel::ComputeDEDXPerVolume ( const G4Material material,
const G4ParticleDefinition ,
G4double  kineticEnergy,
G4double  cutEnergy 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 211 of file G4LivermoreIonisationModel.cc.

References G4VEnergySpectrum::AverageEnergy(), G4VEnergySpectrum::Excitation(), G4VCrossSectionHandler::FindValue(), G4cout, G4endl, G4Material::GetAtomicNumDensityVector(), G4Material::GetElementVector(), G4Material::GetNumberOfElements(), python.hepunit::keV, python.hepunit::mm, n, and G4AtomicTransitionManager::NumberOfShells().

215 {
216  G4double sPower = 0.0;
217 
218  const G4ElementVector* theElementVector = material->GetElementVector();
219  size_t NumberOfElements = material->GetNumberOfElements() ;
220  const G4double* theAtomicNumDensityVector =
221  material->GetAtomicNumDensityVector();
222 
223  // loop for elements in the material
224  for (size_t iel=0; iel<NumberOfElements; iel++ )
225  {
226  G4int iZ = (G4int)((*theElementVector)[iel]->GetZ());
227  G4int nShells = transitionManager->NumberOfShells(iZ);
228  for (G4int n=0; n<nShells; n++)
229  {
230  G4double e = energySpectrum->AverageEnergy(iZ, 0.0,cutEnergy,
231  kineticEnergy, n);
232  G4double cs= crossSectionHandler->FindValue(iZ,kineticEnergy, n);
233  sPower += e * cs * theAtomicNumDensityVector[iel];
234  }
235  G4double esp = energySpectrum->Excitation(iZ,kineticEnergy);
236  sPower += esp * theAtomicNumDensityVector[iel];
237  }
238 
239  if (verboseLevel > 2)
240  {
241  G4cout << "G4LivermoreIonisationModel " << G4endl;
242  G4cout << "Stopping power < " << cutEnergy/keV
243  << " keV at " << kineticEnergy/keV << " keV = "
244  << sPower/(keV/mm) << " keV/mm" << G4endl;
245  }
246 
247  return sPower;
248 }
std::vector< G4Element * > G4ElementVector
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
int G4int
Definition: G4Types.hh:78
G4double FindValue(G4int Z, G4double e) const
virtual G4double AverageEnergy(G4int Z, G4double minKineticEnergy, G4double maxKineticEnergy, G4double kineticEnergy, G4int shell=0, const G4ParticleDefinition *pd=0) const =0
G4GLOB_DLL std::ostream G4cout
const G4int n
const G4double * GetAtomicNumDensityVector() const
Definition: G4Material.hh:214
virtual G4double Excitation(G4int Z, G4double kineticEnergy) const =0
#define G4endl
Definition: G4ios.hh:61
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
double G4double
Definition: G4Types.hh:76
G4int G4LivermoreIonisationModel::GetVerboseLevel ( )
inline

Definition at line 90 of file G4LivermoreIonisationModel.hh.

90 {return verboseLevel;};
void G4LivermoreIonisationModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector cuts 
)
virtual

Implements G4VEmModel.

Definition at line 102 of file G4LivermoreIonisationModel.cc.

References G4VCrossSectionHandler::BuildMeanFreePathForMaterials(), G4VCrossSectionHandler::Clear(), G4Electron::Electron(), FatalException, fParticleChange, G4cout, G4endl, G4Exception(), G4VEmModel::GetParticleChangeForLoss(), python.hepunit::GeV, G4VEmModel::HighEnergyLimit(), python.hepunit::keV, G4VCrossSectionHandler::LoadShellData(), G4VEmModel::LowEnergyLimit(), G4VEnergySpectrum::PrintData(), and G4VCrossSectionHandler::PrintData().

104 {
105  //Check that the Livermore Ionisation is NOT attached to e+
106  if (particle != G4Electron::Electron())
107  {
108  G4Exception("G4LivermoreIonisationModel::Initialise",
109  "em0002",FatalException,
110  "Livermore Ionisation Model is applicable only to electrons");
111  }
112 
113  //Read energy spectrum
114  if (energySpectrum)
115  {
116  delete energySpectrum;
117  energySpectrum = 0;
118  }
119  energySpectrum = new G4eIonisationSpectrum();
120  if (verboseLevel > 3)
121  G4cout << "G4VEnergySpectrum is initialized" << G4endl;
122 
123  //Initialize cross section handler
124  if (crossSectionHandler)
125  {
126  delete crossSectionHandler;
127  crossSectionHandler = 0;
128  }
129 
130  const size_t nbins = 20;
131  G4double emin = LowEnergyLimit();
132  G4double emax = HighEnergyLimit();
133  G4int ndec = G4int(std::log10(emax/emin) + 0.5);
134  if(ndec <= 0) { ndec = 1; }
135 
136  G4VDataSetAlgorithm* interpolation = new G4SemiLogInterpolation();
137  crossSectionHandler =
138  new G4eIonisationCrossSectionHandler(energySpectrum,interpolation,
139  emin,emax,nbins*ndec);
140  crossSectionHandler->Clear();
141  crossSectionHandler->LoadShellData("ioni/ion-ss-cs-");
142  //This is used to retrieve cross section values later on
143  G4VEMDataSet* emdata =
144  crossSectionHandler->BuildMeanFreePathForMaterials(&cuts);
145  //The method BuildMeanFreePathForMaterials() is required here only to force
146  //the building of an internal table: the output pointer can be deleted
147  delete emdata;
148 
149  if (verboseLevel > 0)
150  {
151  G4cout << "Livermore Ionisation model is initialized " << G4endl
152  << "Energy range: "
153  << LowEnergyLimit() / keV << " keV - "
154  << HighEnergyLimit() / GeV << " GeV"
155  << G4endl;
156  }
157 
158  if (verboseLevel > 3)
159  {
160  G4cout << "Cross section data: " << G4endl;
161  crossSectionHandler->PrintData();
162  G4cout << "Parameters: " << G4endl;
163  energySpectrum->PrintData();
164  }
165 
166  if(isInitialised) { return; }
168  isInitialised = true;
169 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:599
virtual void PrintData() const =0
G4ParticleChangeForLoss * GetParticleChangeForLoss()
Definition: G4VEmModel.cc:107
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
int G4int
Definition: G4Types.hh:78
G4ParticleChangeForLoss * fParticleChange
G4GLOB_DLL std::ostream G4cout
G4VEMDataSet * BuildMeanFreePathForMaterials(const G4DataVector *energyCuts=0)
void LoadShellData(const G4String &dataFile)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void G4LivermoreIonisationModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle aDynamicParticle,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 252 of file G4LivermoreIonisationModel.cc.

References G4InuclSpecialFunctions::bindingEnergy(), G4AtomicShell::BindingEnergy(), G4Electron::Electron(), python.hepunit::electron_mass_c2, python.hepunit::eV, fParticleChange, G4cout, G4endl, G4UniformRand, G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), python.hepunit::keV, G4VEnergySpectrum::MaxEnergyOfSecondaries(), G4INCL::Math::min(), G4VParticleChange::ProposeLocalEnergyDeposit(), G4ParticleChangeForLoss::ProposeMomentumDirection(), CLHEP::Hep3Vector::rotateUz(), G4VEnergySpectrum::SampleEnergy(), G4VCrossSectionHandler::SelectRandomAtom(), G4VCrossSectionHandler::SelectRandomShell(), G4DynamicParticle::SetDefinition(), G4DynamicParticle::SetKineticEnergy(), G4DynamicParticle::SetMomentumDirection(), G4ParticleChangeForLoss::SetProposedKineticEnergy(), G4AtomicTransitionManager::Shell(), python.hepunit::twopi, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

258 {
259 
260  G4double kineticEnergy = aDynamicParticle->GetKineticEnergy();
261 
262  if (kineticEnergy <= fIntrinsicLowEnergyLimit)
263  {
266  return;
267  }
268 
269  // Select atom and shell
270  G4int Z = crossSectionHandler->SelectRandomAtom(couple, kineticEnergy);
271  G4int shellIndex = crossSectionHandler->SelectRandomShell(Z, kineticEnergy);
272  const G4AtomicShell* shell = transitionManager->Shell(Z,shellIndex);
274 
275  // Sample delta energy using energy interval for delta-electrons
276  G4double energyMax =
277  std::min(maxE,energySpectrum->MaxEnergyOfSecondaries(kineticEnergy));
278  G4double energyDelta = energySpectrum->SampleEnergy(Z, cutE, energyMax,
279  kineticEnergy, shellIndex);
280 
281  if (energyDelta == 0.) //nothing happens
282  { return; }
283 
284  // Transform to shell potential
285  G4double deltaKinE = energyDelta + 2.0*bindingEnergy;
286  G4double primaryKinE = kineticEnergy + 2.0*bindingEnergy;
287 
288  // sampling of scattering angle neglecting atomic motion
289  G4double deltaMom = std::sqrt(deltaKinE*(deltaKinE + 2.0*electron_mass_c2));
290  G4double primaryMom = std::sqrt(primaryKinE*(primaryKinE + 2.0*electron_mass_c2));
291 
292  G4double cost = deltaKinE * (primaryKinE + 2.0*electron_mass_c2)
293  / (deltaMom * primaryMom);
294  if (cost > 1.) { cost = 1.; }
295  G4double sint = std::sqrt((1. - cost)*(1. + cost));
296  G4double phi = twopi * G4UniformRand();
297  G4double dirx = sint * std::cos(phi);
298  G4double diry = sint * std::sin(phi);
299  G4double dirz = cost;
300 
301  // Rotate to incident electron direction
302  G4ThreeVector primaryDirection = aDynamicParticle->GetMomentumDirection();
303  G4ThreeVector deltaDir(dirx,diry,dirz);
304  deltaDir.rotateUz(primaryDirection);
305  //Updated components
306  dirx = deltaDir.x();
307  diry = deltaDir.y();
308  dirz = deltaDir.z();
309 
310  // Take into account atomic motion del is relative momentum of the motion
311  // kinetic energy of the motion == bindingEnergy in V.Ivanchenko model
312  cost = 2.0*G4UniformRand() - 1.0;
313  sint = std::sqrt(1. - cost*cost);
314  phi = twopi * G4UniformRand();
315  G4double del = std::sqrt(bindingEnergy *(bindingEnergy + 2.0*electron_mass_c2))
316  / deltaMom;
317  dirx += del* sint * std::cos(phi);
318  diry += del* sint * std::sin(phi);
319  dirz += del* cost;
320 
321  // Find out new primary electron direction
322  G4double finalPx = primaryMom*primaryDirection.x() - deltaMom*dirx;
323  G4double finalPy = primaryMom*primaryDirection.y() - deltaMom*diry;
324  G4double finalPz = primaryMom*primaryDirection.z() - deltaMom*dirz;
325 
326  //Ok, ready to create the delta ray
327  G4DynamicParticle* theDeltaRay = new G4DynamicParticle();
328  theDeltaRay->SetKineticEnergy(energyDelta);
329  G4double norm = 1.0/std::sqrt(dirx*dirx + diry*diry + dirz*dirz);
330  dirx *= norm;
331  diry *= norm;
332  dirz *= norm;
333  theDeltaRay->SetMomentumDirection(dirx, diry, dirz);
334  theDeltaRay->SetDefinition(G4Electron::Electron());
335  fvect->push_back(theDeltaRay);
336 
337  //This is the amount of energy available for fluorescence
338  G4double theEnergyDeposit = bindingEnergy;
339 
340  // fill ParticleChange
341  // changed energy and momentum of the actual particle
342  G4double finalKinEnergy = kineticEnergy - energyDelta - theEnergyDeposit;
343  if(finalKinEnergy < 0.0)
344  {
345  theEnergyDeposit += finalKinEnergy;
346  finalKinEnergy = 0.0;
347  }
348  else
349  {
350  G4double normLocal = 1.0/std::sqrt(finalPx*finalPx+finalPy*finalPy+finalPz*finalPz);
351  finalPx *= normLocal;
352  finalPy *= normLocal;
353  finalPz *= normLocal;
354  fParticleChange->ProposeMomentumDirection(finalPx, finalPy, finalPz);
355  }
356  fParticleChange->SetProposedKineticEnergy(finalKinEnergy);
357 
358  if (theEnergyDeposit < 0)
359  {
360  G4cout << "G4LivermoreIonisationModel: Negative energy deposit: "
361  << theEnergyDeposit/eV << " eV" << G4endl;
362  theEnergyDeposit = 0.0;
363  }
364 
365  //Assign local energy deposit
366  fParticleChange->ProposeLocalEnergyDeposit(theEnergyDeposit);
367 
368  if (verboseLevel > 1)
369  {
370  G4cout << "-----------------------------------------------------------" << G4endl;
371  G4cout << "Energy balance from G4LivermoreIonisation" << G4endl;
372  G4cout << "Incoming primary energy: " << kineticEnergy/keV << " keV" << G4endl;
373  G4cout << "-----------------------------------------------------------" << G4endl;
374  G4cout << "Outgoing primary energy: " << finalKinEnergy/keV << " keV" << G4endl;
375  G4cout << "Delta ray " << energyDelta/keV << " keV" << G4endl;
376  G4cout << "Fluorescence: " << (bindingEnergy-theEnergyDeposit)/keV << " keV" << G4endl;
377  G4cout << "Local energy deposit " << theEnergyDeposit/keV << " keV" << G4endl;
378  G4cout << "Total final state: " << (finalKinEnergy+energyDelta+bindingEnergy)
379  << " keV" << G4endl;
380  G4cout << "-----------------------------------------------------------" << G4endl;
381  }
382  return;
383 }
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4double GetKineticEnergy() const
double x() const
G4int SelectRandomShell(G4int Z, G4double e) const
void SetMomentumDirection(const G4ThreeVector &aDirection)
G4double BindingEnergy() const
int G4int
Definition: G4Types.hh:78
double z() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
G4ParticleChangeForLoss * fParticleChange
virtual G4double SampleEnergy(G4int Z, G4double minKineticEnergy, G4double maxKineticEnergy, G4double kineticEnergy, G4int shell=0, const G4ParticleDefinition *pd=0) const =0
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
float electron_mass_c2
Definition: hepunit.py:274
virtual G4double MaxEnergyOfSecondaries(G4double kineticEnergy, G4int Z=0, const G4ParticleDefinition *pd=0) const =0
void SetKineticEnergy(G4double aEnergy)
double y() const
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
G4double bindingEnergy(G4int A, G4int Z)
G4AtomicShell * Shell(G4int Z, size_t shellIndex) const
G4int SelectRandomAtom(const G4MaterialCutsCouple *couple, G4double e) const
void G4LivermoreIonisationModel::SetVerboseLevel ( G4int  vl)
inline

Definition at line 89 of file G4LivermoreIonisationModel.hh.

89 {verboseLevel = vl;};

Field Documentation

G4ParticleChangeForLoss* G4LivermoreIonisationModel::fParticleChange
protected

Definition at line 90 of file G4LivermoreIonisationModel.hh.

Referenced by Initialise(), and SampleSecondaries().


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