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

#include <G4LowEnergyIonisation.hh>

Inheritance diagram for G4LowEnergyIonisation:
G4eLowEnergyLoss G4RDVeLowEnergyLoss G4VContinuousDiscreteProcess G4VProcess

Public Member Functions

 G4LowEnergyIonisation (const G4String &processName="LowEnergyIoni")
 
 ~G4LowEnergyIonisation ()
 
G4bool IsApplicable (const G4ParticleDefinition &)
 
void PrintInfoDefinition ()
 
void BuildPhysicsTable (const G4ParticleDefinition &ParticleType)
 
G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &step)
 
void SetCutForLowEnSecPhotons (G4double cut)
 
void SetCutForLowEnSecElectrons (G4double cut)
 
void ActivateAuger (G4bool val)
 
- Public Member Functions inherited from G4eLowEnergyLoss
 G4eLowEnergyLoss (const G4String &)
 
 ~G4eLowEnergyLoss ()
 
void BuildDEDXTable (const G4ParticleDefinition &aParticleType)
 
G4double GetContinuousStepLimit (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)
 
G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &Step)
 
void ActivateFluorescence (G4bool val)
 
G4bool Fluorescence () const
 
- Public Member Functions inherited from G4RDVeLowEnergyLoss
 G4RDVeLowEnergyLoss (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4RDVeLowEnergyLoss (G4RDVeLowEnergyLoss &)
 
virtual ~G4RDVeLowEnergyLoss ()
 
- Public Member Functions inherited from G4VContinuousDiscreteProcess
 G4VContinuousDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VContinuousDiscreteProcess (G4VContinuousDiscreteProcess &)
 
virtual ~G4VContinuousDiscreteProcess ()
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety, G4GPILSelection *selection)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
- Public Member Functions inherited from G4VProcess
 G4VProcess (const G4String &aName="NoName", G4ProcessType aType=fNotDefined)
 
 G4VProcess (const G4VProcess &right)
 
virtual ~G4VProcess ()
 
G4int operator== (const G4VProcess &right) const
 
G4int operator!= (const G4VProcess &right) const
 
G4double GetCurrentInteractionLength () const
 
void SetPILfactor (G4double value)
 
G4double GetPILfactor () const
 
G4double AlongStepGPIL (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)
 
G4double AtRestGPIL (const G4Track &track, G4ForceCondition *condition)
 
G4double PostStepGPIL (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &, G4bool)
 
const G4StringGetPhysicsTableFileName (const G4ParticleDefinition *, const G4String &directory, const G4String &tableName, G4bool ascii=false)
 
const G4StringGetProcessName () const
 
G4ProcessType GetProcessType () const
 
void SetProcessType (G4ProcessType)
 
G4int GetProcessSubType () const
 
void SetProcessSubType (G4int)
 
virtual void StartTracking (G4Track *)
 
virtual void EndTracking ()
 
virtual void SetProcessManager (const G4ProcessManager *)
 
virtual const G4ProcessManagerGetProcessManager ()
 
virtual void ResetNumberOfInteractionLengthLeft ()
 
G4double GetNumberOfInteractionLengthLeft () const
 
G4double GetTotalNumberOfInteractionLengthTraversed () const
 
G4bool isAtRestDoItIsEnabled () const
 
G4bool isAlongStepDoItIsEnabled () const
 
G4bool isPostStepDoItIsEnabled () const
 
virtual void DumpInfo () const
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 
virtual void SetMasterProcess (G4VProcess *masterP)
 
const G4VProcessGetMasterProcess () const
 
virtual void BuildWorkerPhysicsTable (const G4ParticleDefinition &part)
 
virtual void PrepareWorkerPhysicsTable (const G4ParticleDefinition &)
 

Protected Member Functions

G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual std::vector
< G4DynamicParticle * > * 
DeexciteAtom (const G4MaterialCutsCouple *couple, G4double incidentEnergy, G4double eLoss)
 
- Protected Member Functions inherited from G4RDVeLowEnergyLoss
G4double GetLossWithFluct (const G4DynamicParticle *aParticle, const G4MaterialCutsCouple *couple, G4double MeanLoss, G4double step)
 
- Protected Member Functions inherited from G4VContinuousDiscreteProcess
void SetGPILSelection (G4GPILSelection selection)
 
G4GPILSelection GetGPILSelection () const
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4eLowEnergyLoss
static void SetNbOfProcesses (G4int nb)
 
static void PlusNbOfProcesses ()
 
static void MinusNbOfProcesses ()
 
static G4int GetNbOfProcesses ()
 
static void SetLowerBoundEloss (G4double val)
 
static void SetUpperBoundEloss (G4double val)
 
static void SetNbinEloss (G4int nb)
 
static G4double GetLowerBoundEloss ()
 
static G4double GetUpperBoundEloss ()
 
static G4int GetNbinEloss ()
 
- Static Public Member Functions inherited from G4RDVeLowEnergyLoss
static void SetRndmStep (G4bool value)
 
static void SetEnlossFluc (G4bool value)
 
static void SetStepFunction (G4double c1, G4double c2)
 
- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Static Protected Member Functions inherited from G4RDVeLowEnergyLoss
static G4PhysicsTableBuildRangeTable (G4PhysicsTable *theDEDXTable, G4PhysicsTable *theRangeTable, G4double Tmin, G4double Tmax, G4int nbin)
 
static G4PhysicsTableBuildLabTimeTable (G4PhysicsTable *theDEDXTable, G4PhysicsTable *theLabTimeTable, G4double Tmin, G4double Tmax, G4int nbin)
 
static G4PhysicsTableBuildProperTimeTable (G4PhysicsTable *theDEDXTable, G4PhysicsTable *ProperTimeTable, G4double Tmin, G4double Tmax, G4int nbin)
 
static G4PhysicsTableBuildRangeCoeffATable (G4PhysicsTable *theRangeTable, G4PhysicsTable *theCoeffATable, G4double Tmin, G4double Tmax, G4int nbin)
 
static G4PhysicsTableBuildRangeCoeffBTable (G4PhysicsTable *theRangeTable, G4PhysicsTable *theCoeffBTable, G4double Tmin, G4double Tmax, G4int nbin)
 
static G4PhysicsTableBuildRangeCoeffCTable (G4PhysicsTable *theRangeTable, G4PhysicsTable *theCoeffCTable, G4double Tmin, G4double Tmax, G4int nbin)
 
static G4PhysicsTableBuildInverseRangeTable (G4PhysicsTable *theRangeTable, G4PhysicsTable *theRangeCoeffATable, G4PhysicsTable *theRangeCoeffBTable, G4PhysicsTable *theRangeCoeffCTable, G4PhysicsTable *theInverseRangeTable, G4double Tmin, G4double Tmax, G4int nbin)
 
- Protected Attributes inherited from G4eLowEnergyLoss
G4PhysicsTabletheLossTable
 
G4double MinKineticEnergy
 
G4double Charge
 
G4double lastCharge
 
- Protected Attributes inherited from G4RDVeLowEnergyLoss
const G4MateriallastMaterial
 
G4int imat
 
G4double f1Fluct
 
G4double f2Fluct
 
G4double e1Fluct
 
G4double e2Fluct
 
G4double rateFluct
 
G4double ipotFluct
 
G4double e1LogFluct
 
G4double e2LogFluct
 
G4double ipotLogFluct
 
const G4int nmaxCont1
 
const G4int nmaxCont2
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChangepParticleChange
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft
 
G4double currentInteractionLength
 
G4double theInitialNumberOfInteractionLength
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType
 
G4int theProcessSubType
 
G4double thePILfactor
 
G4bool enableAtRestDoIt
 
G4bool enableAlongStepDoIt
 
G4bool enablePostStepDoIt
 
G4int verboseLevel
 
- Static Protected Attributes inherited from G4eLowEnergyLoss
static G4PhysicsTabletheDEDXElectronTable = 0
 
static G4PhysicsTabletheDEDXPositronTable = 0
 
static G4PhysicsTabletheRangeElectronTable = 0
 
static G4PhysicsTabletheRangePositronTable = 0
 
static G4PhysicsTabletheInverseRangeElectronTable = 0
 
static G4PhysicsTabletheInverseRangePositronTable = 0
 
static G4PhysicsTabletheLabTimeElectronTable = 0
 
static G4PhysicsTabletheLabTimePositronTable = 0
 
static G4PhysicsTabletheProperTimeElectronTable = 0
 
static G4PhysicsTabletheProperTimePositronTable = 0
 
static G4int NbOfProcesses = 2
 
static G4int CounterOfElectronProcess = 0
 
static G4int CounterOfPositronProcess = 0
 
static G4PhysicsTable ** RecorderOfElectronProcess
 
static G4PhysicsTable ** RecorderOfPositronProcess
 
- Static Protected Attributes inherited from G4RDVeLowEnergyLoss
static G4double ParticleMass
 
static G4double taulow
 
static G4double tauhigh
 
static G4double ltaulow
 
static G4double ltauhigh
 
static G4double dRoverRange = 20*perCent
 
static G4double finalRange = 200*micrometer
 
static G4double c1lim = dRoverRange
 
static G4double c2lim = 2.*(1.-dRoverRange)*finalRange
 
static G4double c3lim = -(1.-dRoverRange)*finalRange*finalRange
 
static G4bool rndmStepFlag = false
 
static G4bool EnlossFlucFlag = true
 

Detailed Description

Definition at line 71 of file G4LowEnergyIonisation.hh.

Constructor & Destructor Documentation

G4LowEnergyIonisation::G4LowEnergyIonisation ( const G4String processName = "LowEnergyIoni")

Definition at line 125 of file G4LowEnergyIonisation.cc.

References python.hepunit::eV, and G4VProcess::verboseLevel.

126  : G4eLowEnergyLoss(nam),
127  crossSectionHandler(0),
128  theMeanFreePath(0),
129  energySpectrum(0),
130  shellVacancy(0)
131 {
132  cutForPhotons = 250.0*eV;
133  cutForElectrons = 250.0*eV;
134  verboseLevel = 0;
135 }
G4int verboseLevel
Definition: G4VProcess.hh:368
G4eLowEnergyLoss(const G4String &)
G4LowEnergyIonisation::~G4LowEnergyIonisation ( )

Definition at line 138 of file G4LowEnergyIonisation.cc.

139 {
140  delete crossSectionHandler;
141  delete energySpectrum;
142  delete theMeanFreePath;
143  delete shellVacancy;
144 }

Member Function Documentation

void G4LowEnergyIonisation::ActivateAuger ( G4bool  val)

Definition at line 713 of file G4LowEnergyIonisation.cc.

References G4RDAtomicDeexcitation::ActivateAugerElectronProduction().

714 {
715  deexcitationManager.ActivateAugerElectronProduction(val);
716 }
void ActivateAugerElectronProduction(G4bool val)
void G4LowEnergyIonisation::BuildPhysicsTable ( const G4ParticleDefinition ParticleType)
virtual

Reimplemented from G4VProcess.

Definition at line 147 of file G4LowEnergyIonisation.cc.

References G4eLowEnergyLoss::BuildDEDXTable(), G4eLowEnergyLoss::CounterOfElectronProcess, G4eLowEnergyLoss::CounterOfPositronProcess, G4Electron::Electron(), G4cout, G4endl, G4eLowEnergyLoss::GetLowerBoundEloss(), G4eLowEnergyLoss::GetNbinEloss(), G4VProcess::GetProcessName(), G4eLowEnergyLoss::GetUpperBoundEloss(), G4RDVCrossSectionHandler::LoadShellData(), G4RDVEMDataSet::PrintData(), G4RDVEnergySpectrum::PrintData(), G4RDVCrossSectionHandler::PrintData(), PrintInfoDefinition(), G4eLowEnergyLoss::RecorderOfElectronProcess, G4eLowEnergyLoss::RecorderOfPositronProcess, and G4VProcess::verboseLevel.

148 {
149  if(verboseLevel > 0) {
150  G4cout << "G4LowEnergyIonisation::BuildPhysicsTable start"
151  << G4endl;
152  }
153 
154  cutForDelta.clear();
155 
156  // Create and fill IonisationParameters once
157  if( energySpectrum != 0 ) delete energySpectrum;
158  energySpectrum = new G4RDeIonisationSpectrum();
159 
160  if(verboseLevel > 0) {
161  G4cout << "G4RDVEnergySpectrum is initialized"
162  << G4endl;
163  }
164 
165  // Create and fill G4RDCrossSectionHandler once
166 
167  if ( crossSectionHandler != 0 ) delete crossSectionHandler;
168  G4RDVDataSetAlgorithm* interpolation = new G4RDSemiLogInterpolation();
169  G4double lowKineticEnergy = GetLowerBoundEloss();
170  G4double highKineticEnergy = GetUpperBoundEloss();
171  G4int totBin = GetNbinEloss();
172  crossSectionHandler = new G4RDeIonisationCrossSectionHandler(energySpectrum,
173  interpolation,
174  lowKineticEnergy,
175  highKineticEnergy,
176  totBin);
177  crossSectionHandler->LoadShellData("ioni/ion-ss-cs-");
178 
179  if (verboseLevel > 0) {
181  << " is created; Cross section data: "
182  << G4endl;
183  crossSectionHandler->PrintData();
184  G4cout << "Parameters: "
185  << G4endl;
186  energySpectrum->PrintData();
187  }
188 
189  // Build loss table for IonisationIV
190 
191  BuildLossTable(aParticleType);
192 
193  if(verboseLevel > 0) {
194  G4cout << "The loss table is built"
195  << G4endl;
196  }
197 
198  if (&aParticleType==G4Electron::Electron()) {
199 
200  RecorderOfElectronProcess[CounterOfElectronProcess] = (*this).theLossTable;
203 
204  } else {
205 
206  RecorderOfPositronProcess[CounterOfPositronProcess] = (*this).theLossTable;
208  }
209 
210  // Build mean free path data using cut values
211 
212  if( theMeanFreePath ) delete theMeanFreePath;
213  theMeanFreePath = crossSectionHandler->
214  BuildMeanFreePathForMaterials(&cutForDelta);
215 
216  if(verboseLevel > 0) {
217  G4cout << "The MeanFreePath table is built"
218  << G4endl;
219  if(verboseLevel > 1) theMeanFreePath->PrintData();
220  }
221 
222  // Build common DEDX table for all ionisation processes
223 
224  BuildDEDXTable(aParticleType);
225 
226  if (verboseLevel > 0) {
227  G4cout << "G4LowEnergyIonisation::BuildPhysicsTable end"
228  << G4endl;
229  }
230 }
G4int verboseLevel
Definition: G4VProcess.hh:368
static G4int CounterOfElectronProcess
void LoadShellData(const G4String &dataFile)
static G4double GetLowerBoundEloss()
int G4int
Definition: G4Types.hh:78
G4GLOB_DLL std::ostream G4cout
void BuildDEDXTable(const G4ParticleDefinition &aParticleType)
virtual void PrintData(void) const =0
virtual void PrintData() const =0
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
static G4PhysicsTable ** RecorderOfPositronProcess
static G4double GetUpperBoundEloss()
static G4PhysicsTable ** RecorderOfElectronProcess
static G4int GetNbinEloss()
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
static G4int CounterOfPositronProcess
double G4double
Definition: G4Types.hh:76
std::vector< G4DynamicParticle * > * G4LowEnergyIonisation::DeexciteAtom ( const G4MaterialCutsCouple couple,
G4double  incidentEnergy,
G4double  eLoss 
)
protectedvirtual

Reimplemented from G4eLowEnergyLoss.

Definition at line 606 of file G4LowEnergyIonisation.cc.

References G4RDAtomicShell::BindingEnergy(), G4Electron::Electron(), G4Gamma::Gamma(), G4RDShellVacancy::GenerateNumberOfIonisations(), G4RDAtomicDeexcitation::GenerateParticles(), G4DynamicParticle::GetDefinition(), G4Material::GetElementVector(), G4DynamicParticle::GetKineticEnergy(), G4MaterialCutsCouple::GetMaterial(), G4Material::GetNumberOfElements(), G4RDAtomicTransitionManager::Instance(), n, position, G4RDVCrossSectionHandler::SelectRandomShell(), G4RDAtomicTransitionManager::Shell(), and G4RDAtomicShell::ShellId().

609 {
610  // create vector of secondary particles
611  const G4Material* material = couple->GetMaterial();
612 
613  std::vector<G4DynamicParticle*>* partVector =
614  new std::vector<G4DynamicParticle*>;
615 
616  if(eLoss > cutForPhotons && eLoss > cutForElectrons) {
617 
618  const G4RDAtomicTransitionManager* transitionManager =
620 
621  size_t nElements = material->GetNumberOfElements();
622  const G4ElementVector* theElementVector = material->GetElementVector();
623 
624  std::vector<G4DynamicParticle*>* secVector = 0;
625  G4DynamicParticle* aSecondary = 0;
626  G4ParticleDefinition* type = 0;
627  G4double e;
629  G4int shell, shellId;
630 
631  // sample secondaries
632 
633  G4double eTot = 0.0;
634  std::vector<G4int> n =
635  shellVacancy->GenerateNumberOfIonisations(couple,
636  incidentEnergy,eLoss);
637  for (size_t i=0; i<nElements; i++) {
638 
639  G4int Z = (G4int)((*theElementVector)[i]->GetZ());
640  size_t nVacancies = n[i];
641 
642  G4double maxE = transitionManager->Shell(Z, 0)->BindingEnergy();
643 
644  if (nVacancies && Z > 5 && (maxE>cutForPhotons || maxE>cutForElectrons)) {
645 
646  for (size_t j=0; j<nVacancies; j++) {
647 
648  shell = crossSectionHandler->SelectRandomShell(Z, incidentEnergy);
649  shellId = transitionManager->Shell(Z, shell)->ShellId();
650  G4double maxEShell =
651  transitionManager->Shell(Z, shell)->BindingEnergy();
652 
653  if (maxEShell>cutForPhotons || maxEShell>cutForElectrons ) {
654 
655  secVector = deexcitationManager.GenerateParticles(Z, shellId);
656 
657  if (secVector != 0) {
658 
659  for (size_t l = 0; l<secVector->size(); l++) {
660 
661  aSecondary = (*secVector)[l];
662  if (aSecondary != 0) {
663 
664  e = aSecondary->GetKineticEnergy();
665  type = aSecondary->GetDefinition();
666  if ( eTot + e <= eLoss &&
667  ((type == G4Gamma::Gamma() && e>cutForPhotons ) ||
668  (type == G4Electron::Electron() && e>cutForElectrons))) {
669 
670  eTot += e;
671  partVector->push_back(aSecondary);
672 
673  } else {
674 
675  delete aSecondary;
676 
677  }
678  }
679  }
680  delete secVector;
681  }
682  }
683  }
684  }
685  }
686  }
687  return partVector;
688 }
std::vector< G4DynamicParticle * > * GenerateParticles(G4int Z, G4int shellId)
std::vector< G4Element * > G4ElementVector
G4double GetKineticEnergy() const
G4ParticleDefinition * GetDefinition() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
int G4int
Definition: G4Types.hh:78
static G4RDAtomicTransitionManager * Instance()
std::vector< G4int > GenerateNumberOfIonisations(const G4MaterialCutsCouple *couple, G4double incidentEnergy, G4double eLoss) const
string material
Definition: eplot.py:19
G4double BindingEnergy() const
G4int SelectRandomShell(G4int Z, G4double e) const
const G4int n
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
int position
Definition: filter.cc:7
G4RDAtomicShell * Shell(G4int Z, size_t shellIndex) const
static G4Electron * Electron()
Definition: G4Electron.cc:94
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
double G4double
Definition: G4Types.hh:76
const G4Material * GetMaterial() const
G4int ShellId() const
G4double G4LowEnergyIonisation::GetMeanFreePath ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
protectedvirtual

Implements G4eLowEnergyLoss.

Definition at line 690 of file G4LowEnergyIonisation.cc.

References G4RDVEMDataSet::FindValue(), G4RDVEMDataSet::GetComponent(), G4Track::GetKineticEnergy(), G4Track::GetMaterialCutsCouple(), and NotForced.

693 {
694  *cond = NotForced;
695  G4int index = (track.GetMaterialCutsCouple())->GetIndex();
696  const G4RDVEMDataSet* data = theMeanFreePath->GetComponent(index);
697  G4double meanFreePath = data->FindValue(track.GetKineticEnergy());
698  return meanFreePath;
699 }
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
virtual const G4RDVEMDataSet * GetComponent(G4int componentId) const =0
int G4int
Definition: G4Types.hh:78
G4double GetKineticEnergy() const
virtual G4double FindValue(G4double x, G4int componentId=0) const =0
double G4double
Definition: G4Types.hh:76
const XML_Char const XML_Char * data
G4bool G4LowEnergyIonisation::IsApplicable ( const G4ParticleDefinition particle)
virtual

Reimplemented from G4eLowEnergyLoss.

Definition at line 600 of file G4LowEnergyIonisation.cc.

References G4Electron::Electron().

601 {
602  return ( (&particle == G4Electron::Electron()) );
603 }
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4VParticleChange * G4LowEnergyIonisation::PostStepDoIt ( const G4Track track,
const G4Step step 
)
virtual

Implements G4eLowEnergyLoss.

Definition at line 405 of file G4LowEnergyIonisation.cc.

References G4ParticleChange::AddSecondary(), G4VProcess::aParticleChange, G4InuclSpecialFunctions::bindingEnergy(), G4RDAtomicShell::BindingEnergy(), G4Electron::Electron(), python.hepunit::electron_mass_c2, python.hepunit::eV, G4eLowEnergyLoss::Fluorescence(), fStopAndKill, G4cout, G4endl, G4UniformRand, G4Gamma::Gamma(), G4RDAtomicDeexcitation::GenerateParticles(), G4DynamicParticle::GetDefinition(), G4MaterialCutsCouple::GetIndex(), G4DynamicParticle::GetKineticEnergy(), G4Track::GetKineticEnergy(), G4Track::GetMaterialCutsCouple(), G4Track::GetMomentumDirection(), G4ParticleChange::Initialize(), G4RDAtomicTransitionManager::Instance(), G4RDVEnergySpectrum::MaxEnergyOfSecondaries(), G4VContinuousDiscreteProcess::PostStepDoIt(), G4ParticleChange::ProposeEnergy(), G4VParticleChange::ProposeLocalEnergyDeposit(), G4ParticleChange::ProposeMomentumDirection(), G4VParticleChange::ProposeTrackStatus(), CLHEP::Hep3Vector::rotateUz(), G4RDVEnergySpectrum::SampleEnergy(), G4RDVCrossSectionHandler::SelectRandomAtom(), G4RDVCrossSectionHandler::SelectRandomShell(), G4DynamicParticle::SetDefinition(), G4DynamicParticle::SetKineticEnergy(), G4DynamicParticle::SetMomentumDirection(), G4VParticleChange::SetNumberOfSecondaries(), G4RDAtomicShell::ShellId(), python.hepunit::twopi, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), and CLHEP::Hep3Vector::z().

407 {
408  // Delta electron production mechanism on base of the model
409  // J. Stepanek " A program to determine the radiation spectra due
410  // to a single atomic subshell ionisation by a particle or due to
411  // deexcitation or decay of radionuclides",
412  // Comp. Phys. Comm. 1206 pp 1-19 (1997)
413 
415 
416  const G4MaterialCutsCouple* couple = track.GetMaterialCutsCouple();
417  G4double kineticEnergy = track.GetKineticEnergy();
418 
419  // Select atom and shell
420 
421  G4int Z = crossSectionHandler->SelectRandomAtom(couple, kineticEnergy);
422  G4int shell = crossSectionHandler->SelectRandomShell(Z, kineticEnergy);
423  const G4RDAtomicShell* atomicShell =
424  (G4RDAtomicTransitionManager::Instance())->Shell(Z, shell);
425  G4double bindingEnergy = atomicShell->BindingEnergy();
426  G4int shellId = atomicShell->ShellId();
427 
428  // Sample delta energy
429 
430  G4int index = couple->GetIndex();
431  G4double tCut = cutForDelta[index];
432  G4double tmax = energySpectrum->MaxEnergyOfSecondaries(kineticEnergy);
433  G4double tDelta = energySpectrum->SampleEnergy(Z, tCut, tmax,
434  kineticEnergy, shell);
435 
436  if(tDelta == 0.0)
437  return G4VContinuousDiscreteProcess::PostStepDoIt(track, step);
438 
439  // Transform to shell potential
440  G4double deltaKinE = tDelta + 2.0*bindingEnergy;
441  G4double primaryKinE = kineticEnergy + 2.0*bindingEnergy;
442 
443  // sampling of scattering angle neglecting atomic motion
444  G4double deltaMom = std::sqrt(deltaKinE*(deltaKinE + 2.0*electron_mass_c2));
445  G4double primaryMom = std::sqrt(primaryKinE*(primaryKinE + 2.0*electron_mass_c2));
446 
447  G4double cost = deltaKinE * (primaryKinE + 2.0*electron_mass_c2)
448  / (deltaMom * primaryMom);
449 
450  if (cost > 1.) cost = 1.;
451  G4double sint = std::sqrt(1. - cost*cost);
452  G4double phi = twopi * G4UniformRand();
453  G4double dirx = sint * std::cos(phi);
454  G4double diry = sint * std::sin(phi);
455  G4double dirz = cost;
456 
457  // Rotate to incident electron direction
458  G4ThreeVector primaryDirection = track.GetMomentumDirection();
459  G4ThreeVector deltaDir(dirx,diry,dirz);
460  deltaDir.rotateUz(primaryDirection);
461  dirx = deltaDir.x();
462  diry = deltaDir.y();
463  dirz = deltaDir.z();
464 
465 
466  // Take into account atomic motion del is relative momentum of the motion
467  // kinetic energy of the motion == bindingEnergy in V.Ivanchenko model
468 
469  cost = 2.0*G4UniformRand() - 1.0;
470  sint = std::sqrt(1. - cost*cost);
471  phi = twopi * G4UniformRand();
472  G4double del = std::sqrt(bindingEnergy *(bindingEnergy + 2.0*electron_mass_c2))
473  / deltaMom;
474  dirx += del* sint * std::cos(phi);
475  diry += del* sint * std::sin(phi);
476  dirz += del* cost;
477 
478  // Find out new primary electron direction
479  G4double finalPx = primaryMom*primaryDirection.x() - deltaMom*dirx;
480  G4double finalPy = primaryMom*primaryDirection.y() - deltaMom*diry;
481  G4double finalPz = primaryMom*primaryDirection.z() - deltaMom*dirz;
482 
483  // create G4DynamicParticle object for delta ray
484  G4DynamicParticle* theDeltaRay = new G4DynamicParticle();
485  theDeltaRay->SetKineticEnergy(tDelta);
486  G4double norm = 1.0/std::sqrt(dirx*dirx + diry*diry + dirz*dirz);
487  dirx *= norm;
488  diry *= norm;
489  dirz *= norm;
490  theDeltaRay->SetMomentumDirection(dirx, diry, dirz);
491  theDeltaRay->SetDefinition(G4Electron::Electron());
492 
493  G4double theEnergyDeposit = bindingEnergy;
494 
495  // fill ParticleChange
496  // changed energy and momentum of the actual particle
497 
498  G4double finalKinEnergy = kineticEnergy - tDelta - theEnergyDeposit;
499  if(finalKinEnergy < 0.0) {
500  theEnergyDeposit += finalKinEnergy;
501  finalKinEnergy = 0.0;
503 
504  } else {
505 
506  G4double norm = 1.0/std::sqrt(finalPx*finalPx+finalPy*finalPy+finalPz*finalPz);
507  finalPx *= norm;
508  finalPy *= norm;
509  finalPz *= norm;
510  aParticleChange.ProposeMomentumDirection(finalPx, finalPy, finalPz);
511  }
512 
513  aParticleChange.ProposeEnergy(finalKinEnergy);
514 
515  // Generation of Fluorescence and Auger
516  size_t nSecondaries = 0;
517  size_t totalNumber = 1;
518  std::vector<G4DynamicParticle*>* secondaryVector = 0;
519  G4DynamicParticle* aSecondary = 0;
520  G4ParticleDefinition* type = 0;
521 
522  // Fluorescence data start from element 6
523 
524  if (Fluorescence() && Z > 5 && (bindingEnergy >= cutForPhotons
525  || bindingEnergy >= cutForElectrons)) {
526 
527  secondaryVector = deexcitationManager.GenerateParticles(Z, shellId);
528 
529  if (secondaryVector != 0) {
530 
531  nSecondaries = secondaryVector->size();
532  for (size_t i = 0; i<nSecondaries; i++) {
533 
534  aSecondary = (*secondaryVector)[i];
535  if (aSecondary) {
536 
537  G4double e = aSecondary->GetKineticEnergy();
538  type = aSecondary->GetDefinition();
539  if (e < theEnergyDeposit &&
540  ((type == G4Gamma::Gamma() && e > cutForPhotons ) ||
541  (type == G4Electron::Electron() && e > cutForElectrons ))) {
542 
543  theEnergyDeposit -= e;
544  totalNumber++;
545 
546  } else {
547 
548  delete aSecondary;
549  (*secondaryVector)[i] = 0;
550  }
551  }
552  }
553  }
554  }
555 
556  // Save delta-electrons
557 
559  aParticleChange.AddSecondary(theDeltaRay);
560 
561  // Save Fluorescence and Auger
562 
563  if (secondaryVector) {
564 
565  for (size_t l = 0; l < nSecondaries; l++) {
566 
567  aSecondary = (*secondaryVector)[l];
568 
569  if(aSecondary) {
570 
571  aParticleChange.AddSecondary(aSecondary);
572  }
573  }
574  delete secondaryVector;
575  }
576 
577  if(theEnergyDeposit < 0.) {
578  G4cout << "G4LowEnergyIonisation: Negative energy deposit: "
579  << theEnergyDeposit/eV << " eV" << G4endl;
580  theEnergyDeposit = 0.0;
581  }
582  aParticleChange.ProposeLocalEnergyDeposit(theEnergyDeposit);
583 
584  return G4VContinuousDiscreteProcess::PostStepDoIt(track, step);
585 }
std::vector< G4DynamicParticle * > * GenerateParticles(G4int Z, G4int shellId)
G4double GetKineticEnergy() const
double x() const
G4bool Fluorescence() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
void SetMomentumDirection(const G4ThreeVector &aDirection)
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
double z() const
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
static G4RDAtomicTransitionManager * Instance()
G4double GetKineticEnergy() const
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
G4double BindingEnergy() const
G4int SelectRandomShell(G4int Z, G4double e) const
float electron_mass_c2
Definition: hepunit.py:274
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void SetKineticEnergy(G4double aEnergy)
virtual void Initialize(const G4Track &)
const G4ThreeVector & GetMomentumDirection() const
void SetNumberOfSecondaries(G4int totSecondaries)
double y() const
G4int SelectRandomAtom(const G4MaterialCutsCouple *couple, G4double e) const
void ProposeEnergy(G4double finalEnergy)
G4ParticleChange aParticleChange
Definition: G4VProcess.hh:289
void AddSecondary(G4Track *aSecondary)
static G4Electron * Electron()
Definition: G4Electron.cc:94
#define G4endl
Definition: G4ios.hh:61
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
void SetDefinition(const G4ParticleDefinition *aParticleDefinition)
G4double bindingEnergy(G4int A, G4int Z)
virtual G4double SampleEnergy(G4int Z, G4double minKineticEnergy, G4double maxKineticEnergy, G4double kineticEnergy, G4int shell=0, const G4ParticleDefinition *pd=0) const =0
G4int ShellId() const
virtual G4double MaxEnergyOfSecondaries(G4double kineticEnergy, G4int Z=0, const G4ParticleDefinition *pd=0) const =0
void G4LowEnergyIonisation::PrintInfoDefinition ( )

Definition at line 588 of file G4LowEnergyIonisation.cc.

References G4cout, G4endl, and G4VProcess::GetProcessName().

Referenced by BuildPhysicsTable().

589 {
590  G4String comments = "Total cross sections from EEDL database.";
591  comments += "\n Gamma energy sampled from a parametrised formula.";
592  comments += "\n Implementation of the continuous dE/dx part.";
593  comments += "\n At present it can be used for electrons ";
594  comments += "in the energy range [250eV,100GeV].";
595  comments += "\n The process must work with G4LowEnergyBremsstrahlung.";
596 
597  G4cout << G4endl << GetProcessName() << ": " << comments << G4endl;
598 }
G4GLOB_DLL std::ostream G4cout
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
#define G4endl
Definition: G4ios.hh:61
void G4LowEnergyIonisation::SetCutForLowEnSecElectrons ( G4double  cut)

Definition at line 707 of file G4LowEnergyIonisation.cc.

References G4RDAtomicDeexcitation::SetCutForAugerElectrons().

708 {
709  cutForElectrons = cut;
710  deexcitationManager.SetCutForAugerElectrons(cut);
711 }
void SetCutForAugerElectrons(G4double cut)
void G4LowEnergyIonisation::SetCutForLowEnSecPhotons ( G4double  cut)

Definition at line 701 of file G4LowEnergyIonisation.cc.

References G4RDAtomicDeexcitation::SetCutForSecondaryPhotons().

702 {
703  cutForPhotons = cut;
704  deexcitationManager.SetCutForSecondaryPhotons(cut);
705 }
void SetCutForSecondaryPhotons(G4double cut)

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