Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4ScreenedNuclearRecoil Class Reference

A process which handles screened Coulomb collisions between nuclei. More...

#include <G4ScreenedNuclearRecoil.hh>

Inheritance diagram for G4ScreenedNuclearRecoil:
Collaboration diagram for G4ScreenedNuclearRecoil:

Public Member Functions

 G4ScreenedNuclearRecoil (const G4String &processName="ScreenedElastic", const G4String &ScreeningKey="zbl", G4bool GenerateRecoils=1, G4double RecoilCutoff=100.0 *eV, G4double PhysicsCutoff=10.0 *eV)
 Construct the process and set some physics parameters for it. More...
 
virtual ~G4ScreenedNuclearRecoil ()
 destructor More...
 
virtual G4double GetMeanFreePath (const G4Track &, G4double, G4ForceCondition *)
 used internally by Geant4 machinery More...
 
virtual G4VParticleChangePostStepDoIt (const G4Track &aTrack, const G4Step &aStep)
 used internally by Geant4 machinery More...
 
virtual G4bool IsApplicable (const G4ParticleDefinition &aParticleType)
 test if a prticle of type aParticleType can use this process More...
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &aParticleType)
 Build physics tables in advance. Not Implemented. More...
 
virtual void DumpPhysicsTable (const G4ParticleDefinition &aParticleType)
 Export physics tables for persistency. Not Implemented. More...
 
virtual G4bool CheckNuclearCollision (G4double A, G4double A1, G4double apsis)
 deterine if the moving particle is within the strong force range of the selected nucleus More...
 
virtual
G4ScreenedCoulombCrossSection
GetNewCrossSectionHandler (void)
 
G4double GetNIEL () const
 Get non-ionizing energy loss for last step. More...
 
void ResetTables ()
 clear precomputed screening tables More...
 
void SetMaxEnergyForScattering (G4double energy)
 set the upper energy beyond which this process has no cross section More...
 
std::string GetScreeningKey () const
 find out what screening function we are using More...
 
void AllowEnergyDeposition (G4bool flag)
 enable or disable all energy deposition by this process More...
 
G4bool GetAllowEnergyDeposition () const
 get flag indicating whether deposition is enabled More...
 
void EnableRecoils (G4bool flag)
 enable or disable the generation of recoils. If recoils are disabled, the energy they would have received is just deposited. param flag if true, create recoil ions in cases in which the energy is above the recoilCutoff. If false, just deposit the energy. More...
 
G4bool GetEnableRecoils () const
 find out if generation of recoils is enabled. More...
 
void SetMFPScaling (G4double scale)
 set the mean free path scaling as specified More...
 
G4double GetMFPScaling () const
 get the MFPScaling parameter More...
 
void AvoidNuclearReactions (G4bool flag)
 enable or disable whether this process will skip collisions which are close enough they need hadronic phsyics. Default is true (skip close collisions). Disabling this results in excess nuclear stopping power. More...
 
G4bool GetAvoidNuclearReactions () const
 get the flag indicating whether hadronic collisions are ignored. More...
 
void SetRecoilCutoff (G4double energy)
 set the minimum energy (per nucleon) at which recoils can be generated, and the energy (per nucleon) below which all ions are stopped. More...
 
G4double GetRecoilCutoff () const
 get the recoil cutoff More...
 
void SetPhysicsCutoff (G4double energy)
 set the energy to which screening tables are computed. Typically, this is 10 eV or so, and not often changed. More...
 
G4double GetPhysicsCutoff () const
 get the physics cutoff energy. More...
 
void SetNIELPartitionFunction (const G4VNIELPartition *part)
 set the pointer to a class for paritioning energy into NIEL More...
 
void SetCrossSectionHardening (G4double fraction, G4double HardeningFactor)
 set the cross section boost to provide faster computation of backscattering More...
 
G4double GetHardeningFraction () const
 get the fraction of particles which will have boosted scattering More...
 
G4double GetHardeningFactor () const
 get the boost factor in use. More...
 
G4double GetCurrentInteractionLength () const
 the the interaciton length used in the last scattering. More...
 
void SetExternalCrossSectionHandler (G4ScreenedCoulombCrossSection *cs)
 set a function to compute screening tables, if the user needs non-standard behavior. More...
 
G4int GetVerboseLevel () const
 get the verbosity. More...
 
std::map< G4int,
G4ScreenedCoulombCrossSection * > & 
GetCrossSectionHandlers ()
 
void ClearStages (void)
 
void AddStage (G4ScreenedCollisionStage *stage)
 
G4CoulombKinematicsInfoGetKinematics ()
 
void SetValidCollision (G4bool flag)
 
G4bool GetValidCollision () const
 
class G4ParticleChangeGetParticleChange ()
 get the pointer to our ParticleChange object. for internal use, primarily. More...
 
void DepositEnergy (G4int z1, G4double a1, const G4Material *material, G4double energy)
 take the given energy, and use the material information to partition it into NIEL and ionizing energy. More...
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &, G4ForceCondition *)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &, const G4Step &)
 
virtual G4VParticleChangeAlongStepDoIt (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 Attributes

G4double highEnergyLimit
 the energy per nucleon above which the MFP is constant More...
 
G4double lowEnergyLimit
 the energy per nucleon below which the MFP is zero More...
 
G4double processMaxEnergy
 the energy per nucleon beyond which the cross section is zero, to cross over to G4MSC More...
 
G4String screeningKey
 
G4bool generateRecoils
 
G4bool avoidReactions
 
G4double recoilCutoff
 
G4double physicsCutoff
 
G4bool registerDepositedEnergy
 
G4double IonizingLoss
 
G4double NIEL
 
G4double MFPScale
 
G4double hardeningFraction
 
G4double hardeningFactor
 
G4ScreenedCoulombCrossSectionexternalCrossSectionConstructor
 
std::vector
< G4ScreenedCollisionStage * > 
collisionStages
 
std::map< G4int,
G4ScreenedCoulombCrossSection * > 
crossSectionHandlers
 
G4bool validCollision
 
G4CoulombKinematicsInfo kinematics
 
const G4VNIELPartitionNIELPartitionFunction
 
- 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
 

Friends

class G4ScreenedCollisionStage
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Detailed Description

A process which handles screened Coulomb collisions between nuclei.

Definition at line 208 of file G4ScreenedNuclearRecoil.hh.

Constructor & Destructor Documentation

G4ScreenedNuclearRecoil::G4ScreenedNuclearRecoil ( const G4String processName = "ScreenedElastic",
const G4String ScreeningKey = "zbl",
G4bool  GenerateRecoils = 1,
G4double  RecoilCutoff = 100.0*eV,
G4double  PhysicsCutoff = 10.0*eV 
)

Construct the process and set some physics parameters for it.

Parameters
processNamethe name to assign the process
ScreeningKeythe name of a screening function to use. The default functions are "zbl" (recommended for soft scattering), "lj" (recommended for backscattering) and "mol" (Moliere potential)
GenerateRecoilsif frue, ions struck by primary are converted into new moving particles. If false, energy is deposited, but no new moving ions are created.
RecoilCutoffenergy below which no new moving particles will be created, even if a GenerateRecoils is true. Also, a moving primary particle will be stopped if its energy falls below this limit.
PhysicsCutoffthe energy transfer to which screening tables are calucalted. There is no really compelling reason to change it from the 10.0 eV default. However, see the paper on running this in thin targets for further discussion, and its interaction with SetMFPScaling()

Definition at line 320 of file G4ScreenedNuclearRecoil.cc.

323  :
324  G4VDiscreteProcess(processName, fElectromagnetic),
325  screeningKey(ScreeningKey),
326  generateRecoils(GenerateRecoils), avoidReactions(1),
327  recoilCutoff(RecoilCutoff), physicsCutoff(PhysicsCutoff),
331 {
332  // for now, point to class instance of this. Doing it by creating a new
333  // one fails
334  // to correctly update NIEL
335  // not even this is needed... done in G4VProcess().
336  // pParticleChange=&aParticleChange;
337  processMaxEnergy=50000.0*MeV;
338  highEnergyLimit=100.0*MeV;
340  registerDepositedEnergy=1; // by default, don't hide NIEL
341  MFPScale=1.0;
342  // SetVerboseLevel(2);
346 }
void AddStage(G4ScreenedCollisionStage *stage)
G4double highEnergyLimit
the energy per nucleon above which the MFP is constant
const G4VNIELPartition * NIELPartitionFunction
G4ScreenedCoulombCrossSection * externalCrossSectionConstructor
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
G4double lowEnergyLimit
the energy per nucleon below which the MFP is zero
static constexpr double MeV
Definition: G4SIunits.hh:214
G4double processMaxEnergy
the energy per nucleon beyond which the cross section is zero, to cross over to G4MSC ...

Here is the call graph for this function:

G4ScreenedNuclearRecoil::~G4ScreenedNuclearRecoil ( )
virtual

destructor

Definition at line 391 of file G4ScreenedNuclearRecoil.cc.

392 {
393  ResetTables();
394 }
void ResetTables()
clear precomputed screening tables

Here is the call graph for this function:

Member Function Documentation

void G4ScreenedNuclearRecoil::AddStage ( G4ScreenedCollisionStage stage)
inline

Definition at line 399 of file G4ScreenedNuclearRecoil.hh.

399  {
400  collisionStages.push_back(stage); }
std::vector< G4ScreenedCollisionStage * > collisionStages

Here is the caller graph for this function:

void G4ScreenedNuclearRecoil::AllowEnergyDeposition ( G4bool  flag)
inline

enable or disable all energy deposition by this process

Parameters
flagif true, enable deposition of energy (the default). If false, disable deposition.

Definition at line 296 of file G4ScreenedNuclearRecoil.hh.

void G4ScreenedNuclearRecoil::AvoidNuclearReactions ( G4bool  flag)
inline

enable or disable whether this process will skip collisions which are close enough they need hadronic phsyics. Default is true (skip close collisions). Disabling this results in excess nuclear stopping power.

Parameters
flagtrue results in hard collisions being skipped. false allows hard collisions.

Definition at line 331 of file G4ScreenedNuclearRecoil.hh.

void G4ScreenedNuclearRecoil::BuildPhysicsTable ( const G4ParticleDefinition aParticleType)
virtual

Build physics tables in advance. Not Implemented.

Parameters
aParticleTypethe type of particle to build tables for

Reimplemented from G4VProcess.

Definition at line 805 of file G4ScreenedNuclearRecoil.cc.

806 {
807  G4String nam = aParticleType.GetParticleName();
808  if(nam == "GenericIon" || nam == "proton"
809  || nam == "deuteron" || nam == "triton"
810  || nam == "alpha" || nam == "He3") {
811  G4cout << G4endl << GetProcessName() << ": for " << nam
812  << " SubType= " << GetProcessSubType()
813  << " maxEnergy(MeV)= " << processMaxEnergy/MeV << G4endl;
814  }
815 }
const G4String & GetParticleName() const
G4GLOB_DLL std::ostream G4cout
const G4String & GetProcessName() const
Definition: G4VProcess.hh:408
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
G4int GetProcessSubType() const
Definition: G4VProcess.hh:426
G4double processMaxEnergy
the energy per nucleon beyond which the cross section is zero, to cross over to G4MSC ...

Here is the call graph for this function:

G4bool G4ScreenedNuclearRecoil::CheckNuclearCollision ( G4double  A,
G4double  A1,
G4double  apsis 
)
virtual

deterine if the moving particle is within the strong force range of the selected nucleus

Parameters
Athe nucleon number of the beam
A1the nucleon number of the target
apsisthe distance of closest approach

Definition at line 398 of file G4ScreenedNuclearRecoil.cc.

399  {
400  return avoidReactions && (apsis < (1.1*(std::pow(A,1.0/3.0)+
401  std::pow(a1,1.0/3.0)) + 1.4)*fermi);
402  // nuclei are within 1.4 fm (reduced pion Compton wavelength) of each
403  // other at apsis,
404  // this is hadronic, skip it
405 }
double A(double temperature)
static constexpr double fermi
Definition: G4SIunits.hh:103

Here is the caller graph for this function:

void G4ScreenedNuclearRecoil::ClearStages ( void  )

Definition at line 359 of file G4ScreenedNuclearRecoil.cc.

360 {
361  // I don't think I like deleting the processes here... they are better
362  // abandoned
363  // if the creator doesn't get rid of them
364  // std::vector<G4ScreenedCollisionStage *>::iterator stage=
365  //collisionStages.begin();
366  //for(; stage != collisionStages.end(); stage++) delete (*stage);
367 
368  collisionStages.clear();
369 }
std::vector< G4ScreenedCollisionStage * > collisionStages
void G4ScreenedNuclearRecoil::DepositEnergy ( G4int  z1,
G4double  a1,
const G4Material material,
G4double  energy 
)

take the given energy, and use the material information to partition it into NIEL and ionizing energy.

Definition at line 378 of file G4ScreenedNuclearRecoil.cc.

380 {
381  if(!NIELPartitionFunction) {
383  } else {
385  material, energy);
386  IonizingLoss+=energy*(1-part);
387  NIEL += energy*part;
388  }
389 }
const G4VNIELPartition * NIELPartitionFunction
G4double energy(const ThreeVector &p, const G4double m)
double G4double
Definition: G4Types.hh:76
virtual G4double PartitionNIEL(G4int z1, G4double a1, const G4Material *material, G4double energy) const =0

Here is the call graph for this function:

Here is the caller graph for this function:

void G4ScreenedNuclearRecoil::DumpPhysicsTable ( const G4ParticleDefinition aParticleType)
virtual

Export physics tables for persistency. Not Implemented.

Parameters
aParticleTypethe type of particle to build tables for

Definition at line 819 of file G4ScreenedNuclearRecoil.cc.

820 {
821 }
void G4ScreenedNuclearRecoil::EnableRecoils ( G4bool  flag)
inline

enable or disable the generation of recoils. If recoils are disabled, the energy they would have received is just deposited. param flag if true, create recoil ions in cases in which the energy is above the recoilCutoff. If false, just deposit the energy.

Definition at line 308 of file G4ScreenedNuclearRecoil.hh.

G4bool G4ScreenedNuclearRecoil::GetAllowEnergyDeposition ( ) const
inline

get flag indicating whether deposition is enabled

Definition at line 299 of file G4ScreenedNuclearRecoil.hh.

G4bool G4ScreenedNuclearRecoil::GetAvoidNuclearReactions ( ) const
inline

get the flag indicating whether hadronic collisions are ignored.

Definition at line 334 of file G4ScreenedNuclearRecoil.hh.

std::map<G4int, G4ScreenedCoulombCrossSection*>& G4ScreenedNuclearRecoil::GetCrossSectionHandlers ( )
inline

Definition at line 394 of file G4ScreenedNuclearRecoil.hh.

395  { return crossSectionHandlers; }
std::map< G4int, G4ScreenedCoulombCrossSection * > crossSectionHandlers

Here is the caller graph for this function:

G4double G4ScreenedNuclearRecoil::GetCurrentInteractionLength ( ) const
inline

the the interaciton length used in the last scattering.

Definition at line 380 of file G4ScreenedNuclearRecoil.hh.

380  {
381  return currentInteractionLength; }
G4double currentInteractionLength
Definition: G4VProcess.hh:297

Here is the caller graph for this function:

G4bool G4ScreenedNuclearRecoil::GetEnableRecoils ( ) const
inline

find out if generation of recoils is enabled.

Definition at line 311 of file G4ScreenedNuclearRecoil.hh.

Here is the caller graph for this function:

G4double G4ScreenedNuclearRecoil::GetHardeningFactor ( ) const
inline

get the boost factor in use.

Definition at line 377 of file G4ScreenedNuclearRecoil.hh.

Here is the caller graph for this function:

G4double G4ScreenedNuclearRecoil::GetHardeningFraction ( ) const
inline

get the fraction of particles which will have boosted scattering

Definition at line 374 of file G4ScreenedNuclearRecoil.hh.

Here is the caller graph for this function:

G4CoulombKinematicsInfo& G4ScreenedNuclearRecoil::GetKinematics ( )
inline

Definition at line 402 of file G4ScreenedNuclearRecoil.hh.

402 { return kinematics; }
G4CoulombKinematicsInfo kinematics

Here is the caller graph for this function:

G4double G4ScreenedNuclearRecoil::GetMeanFreePath ( const G4Track track,
G4double  ,
G4ForceCondition cond 
)
virtual

used internally by Geant4 machinery

Implements G4VDiscreteProcess.

Definition at line 417 of file G4ScreenedNuclearRecoil.cc.

420 {
421  const G4DynamicParticle* incoming = track.GetDynamicParticle();
422  G4double energy = incoming->GetKineticEnergy();
423  G4double a1=incoming->GetDefinition()->GetPDGMass()/amu_c2;
424 
425  G4double meanFreePath;
426  *cond=NotForced;
427 
428  if (energy < lowEnergyLimit || energy < recoilCutoff*a1) {
429  *cond=Forced;
430  return 1.0*nm;
431 /* catch and stop slow particles to collect their NIEL! */
432  } else if (energy > processMaxEnergy*a1) {
433  return DBL_MAX; // infinite mean free path
434  } else if (energy > highEnergyLimit*a1) energy=highEnergyLimit*a1;
435 /* constant MFP at high energy */
436 
437  G4double fz1=incoming->GetDefinition()->GetPDGCharge();
438  G4int z1=(G4int)(fz1/eplus + 0.5);
439 
440  std::map<G4int, G4ScreenedCoulombCrossSection*>::iterator xh=
441  crossSectionHandlers.find(z1);
443 
444  if (xh==crossSectionHandlers.end()) {
446  xs->LoadData(screeningKey, z1, a1, physicsCutoff);
447  xs->BuildMFPTables();
448  } else xs=(*xh).second;
449 
450  const G4MaterialCutsCouple* materialCouple =
451  track.GetMaterialCutsCouple();
452  size_t materialIndex = materialCouple->GetMaterial()->GetIndex();
453 
454  const G4_c2_function &mfp=*(*xs)[materialIndex];
455 
456  // make absolutely certain we don't get an out-of-range energy
457  meanFreePath = mfp(std::min(std::max(energy, mfp.xmin()), mfp.xmax()));
458 
459  // G4cout << "MFP: " << meanFreePath << " index " << materialIndex
460  //<< " energy " << energy << " MFPScale " << MFPScale << G4endl;
461 
462  return meanFreePath*MFPScale;
463 }
virtual void LoadData(G4String screeningKey, G4int z1, G4double m1, G4double recoilCutoff)=0
G4double GetKineticEnergy() const
const G4DynamicParticle * GetDynamicParticle() const
G4double highEnergyLimit
the energy per nucleon above which the MFP is constant
size_t GetIndex() const
Definition: G4Material.hh:262
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
std::map< G4int, G4ScreenedCoulombCrossSection * > crossSectionHandlers
static constexpr double eplus
Definition: G4SIunits.hh:199
G4double lowEnergyLimit
the energy per nucleon below which the MFP is zero
G4double GetPDGMass() const
T max(const T t1, const T t2)
brief Return the largest of the two arguments
static constexpr double nm
Definition: G4SIunits.hh:112
G4double energy(const ThreeVector &p, const G4double m)
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
virtual G4ScreenedCoulombCrossSection * GetNewCrossSectionHandler(void)
double G4double
Definition: G4Types.hh:76
float amu_c2
Definition: hepunit.py:277
G4double GetPDGCharge() const
float_type xmax() const
Definition: c2_function.hh:390
#define DBL_MAX
Definition: templates.hh:83
float_type xmin() const
Definition: c2_function.hh:389
const G4Material * GetMaterial() const
G4double processMaxEnergy
the energy per nucleon beyond which the cross section is zero, to cross over to G4MSC ...

Here is the call graph for this function:

G4double G4ScreenedNuclearRecoil::GetMFPScaling ( ) const
inline

get the MFPScaling parameter

Definition at line 322 of file G4ScreenedNuclearRecoil.hh.

322 { return MFPScale; }
G4ScreenedCoulombCrossSection * G4ScreenedNuclearRecoil::GetNewCrossSectionHandler ( void  )
virtual

Definition at line 408 of file G4ScreenedNuclearRecoil.cc.

Here is the call graph for this function:

Here is the caller graph for this function:

G4double G4ScreenedNuclearRecoil::GetNIEL ( ) const
inline

Get non-ionizing energy loss for last step.

Definition at line 270 of file G4ScreenedNuclearRecoil.hh.

270 { return NIEL; }
class G4ParticleChange& G4ScreenedNuclearRecoil::GetParticleChange ( )
inline

get the pointer to our ParticleChange object. for internal use, primarily.

Definition at line 410 of file G4ScreenedNuclearRecoil.hh.

411  { return static_cast<G4ParticleChange &>(*pParticleChange); }

Here is the caller graph for this function:

G4double G4ScreenedNuclearRecoil::GetPhysicsCutoff ( ) const
inline

get the physics cutoff energy.

Definition at line 354 of file G4ScreenedNuclearRecoil.hh.

354 { return physicsCutoff; }
G4double G4ScreenedNuclearRecoil::GetRecoilCutoff ( ) const
inline

get the recoil cutoff

Definition at line 344 of file G4ScreenedNuclearRecoil.hh.

344 { return recoilCutoff; }

Here is the caller graph for this function:

std::string G4ScreenedNuclearRecoil::GetScreeningKey ( ) const
inline

find out what screening function we are using

Definition at line 290 of file G4ScreenedNuclearRecoil.hh.

290 { return screeningKey; }
G4bool G4ScreenedNuclearRecoil::GetValidCollision ( ) const
inline

Definition at line 406 of file G4ScreenedNuclearRecoil.hh.

Here is the caller graph for this function:

G4int G4ScreenedNuclearRecoil::GetVerboseLevel ( ) const
inline

get the verbosity.

Definition at line 392 of file G4ScreenedNuclearRecoil.hh.

392 { return verboseLevel; }
G4int verboseLevel
Definition: G4VProcess.hh:368

Here is the caller graph for this function:

G4bool G4ScreenedNuclearRecoil::IsApplicable ( const G4ParticleDefinition aParticleType)
virtual

test if a prticle of type aParticleType can use this process

Parameters
aParticleTypethe particle to test

Reimplemented from G4VProcess.

Definition at line 796 of file G4ScreenedNuclearRecoil.cc.

797 {
798  return aParticleType == *(G4Proton::Proton()) ||
799  aParticleType.GetParticleType() == "nucleus" ||
800  aParticleType.GetParticleType() == "static_nucleus";
801 }
static G4Proton * Proton()
Definition: G4Proton.cc:93
const G4String & GetParticleType() const

Here is the call graph for this function:

G4VParticleChange * G4ScreenedNuclearRecoil::PostStepDoIt ( const G4Track aTrack,
const G4Step aStep 
)
virtual

used internally by Geant4 machinery

Reimplemented from G4VDiscreteProcess.

Definition at line 465 of file G4ScreenedNuclearRecoil.cc.

467 {
468  validCollision=1;
469  pParticleChange->Initialize(aTrack);
470  NIEL=0.0; // default is no NIEL deposited
471  IonizingLoss=0.0;
472 
473  // do universal setup
474 
475  const G4DynamicParticle* incidentParticle = aTrack.GetDynamicParticle();
476  G4ParticleDefinition *baseParticle=aTrack.GetDefinition();
477 
478  G4double fz1=baseParticle->GetPDGCharge()/eplus;
479  G4int z1=(G4int)(fz1+0.5);
480  G4double a1=baseParticle->GetPDGMass()/amu_c2;
481  G4double incidentEnergy = incidentParticle->GetKineticEnergy();
482 
483  // Select randomly one element and (possibly) isotope in the
484  // current material.
485  const G4MaterialCutsCouple* couple = aTrack.GetMaterialCutsCouple();
486 
487  const G4Material* mat = couple->GetMaterial();
488 
489  G4double P=0.0; // the impact parameter of this collision
490 
491  if(incidentEnergy < GetRecoilCutoff()*a1) {
492  // check energy sanity on entry
493  DepositEnergy(z1, baseParticle->GetPDGMass()/amu_c2, mat,
494  incidentEnergy);
496  // stop the particle and bail out
497  validCollision=0;
498  } else {
499 
500  G4double numberDensity=mat->GetTotNbOfAtomsPerVolume();
501  G4double lattice=0.5/std::pow(numberDensity,1.0/3.0);
502  // typical lattice half-spacing
504  G4double sigopi=1.0/(pi*numberDensity*length);
505  // this is sigma0/pi
506 
507  // compute the impact parameter very early, so if is rejected
508  // as too far away, little effort is wasted
509  // this is the TRIM method for determining an impact parameter
510  // based on the flight path
511  // this gives a cumulative distribution of
512  // N(P)= 1-exp(-pi P^2 n l)
513  // which says the probability of NOT hitting a disk of area
514  // sigma= pi P^2 =exp(-sigma N l)
515  // which may be reasonable
516  if(sigopi < lattice*lattice) {
517  // normal long-flight approximation
518  P = std::sqrt(-std::log(G4UniformRand()) *sigopi);
519  } else {
520  // short-flight limit
521  P = std::sqrt(G4UniformRand())*lattice;
522  }
523 
524  G4double fraction=GetHardeningFraction();
525  if(fraction && G4UniformRand() < fraction) {
526  // pick out some events, and increase the central cross
527  // section by reducing the impact parameter
528  P /= std::sqrt(GetHardeningFactor());
529  }
530 
531 
532  // check if we are far enough away that the energy transfer
533  // must be below cutoff,
534  // and leave everything alone if so, saving a lot of time.
535  if(P*P > sigopi) {
536  if(GetVerboseLevel() > 1)
537  printf("ScreenedNuclear impact reject: length=%.3f P=%.4f limit=%.4f\n",
538  length/angstrom, P/angstrom,std::sqrt(sigopi)/angstrom);
539  // no collision, don't follow up with anything
540  validCollision=0;
541  }
542  }
543 
544  // find out what we hit, and record it in our kinematics block.
546  kinematics.a1=a1;
547 
548  if(validCollision) {
551  G4ParticleDefinition *recoilIon=
552  xsect->SelectRandomUnweightedTarget(couple);
553  kinematics.crossSection=xsect;
554  kinematics.recoilIon=recoilIon;
556  kinematics.a2=recoilIon->GetPDGMass()/amu_c2;
557  } else {
560  kinematics.a2=0;
561  }
562 
563  std::vector<G4ScreenedCollisionStage *>::iterator stage=
564  collisionStages.begin();
565 
566  for(; stage != collisionStages.end(); stage++)
567  (*stage)->DoCollisionStep(this,aTrack, aStep);
568 
572  //MHM G4cout << "depositing energy, total = "
573  //<< IonizingLoss+NIEL << " NIEL = " << NIEL << G4endl;
574  }
575 
576  return G4VDiscreteProcess::PostStepDoIt( aTrack, aStep );
577 }
G4ParticleDefinition * GetDefinition() const
virtual void Initialize(const G4Track &)
G4CoulombKinematicsInfo kinematics
G4double GetCurrentInteractionLength() const
the the interaciton length used in the last scattering.
G4double GetKineticEnergy() const
const G4DynamicParticle * GetDynamicParticle() const
G4ParticleDefinition * recoilIon
G4double GetRecoilCutoff() const
get the recoil cutoff
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
void DepositEnergy(G4int z1, G4double a1, const G4Material *material, G4double energy)
take the given energy, and use the material information to partition it into NIEL and ionizing energy...
int G4int
Definition: G4Types.hh:78
static double P[]
std::vector< G4ScreenedCollisionStage * > collisionStages
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
#define G4UniformRand()
Definition: Randomize.hh:97
G4ScreenedCoulombCrossSection * crossSection
void ProposeNonIonizingEnergyDeposit(G4double anEnergyPart)
G4double GetHardeningFactor() const
get the boost factor in use.
static constexpr double eplus
Definition: G4SIunits.hh:199
G4ParticleDefinition * SelectRandomUnweightedTarget(const G4MaterialCutsCouple *couple)
G4double GetPDGMass() const
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
void ProposeEnergy(G4double finalEnergy)
G4double GetHardeningFraction() const
get the fraction of particles which will have boosted scattering
static constexpr double angstrom
Definition: G4SIunits.hh:102
static constexpr double pi
Definition: G4SIunits.hh:75
double G4double
Definition: G4Types.hh:76
G4int GetVerboseLevel() const
get the verbosity.
float amu_c2
Definition: hepunit.py:277
class G4ParticleChange & GetParticleChange()
get the pointer to our ParticleChange object. for internal use, primarily.
G4double GetPDGCharge() const
virtual G4VParticleChange * PostStepDoIt(const G4Track &, const G4Step &)
const G4Material * GetMaterial() const
std::map< G4int, G4ScreenedCoulombCrossSection * > & GetCrossSectionHandlers()

Here is the call graph for this function:

void G4ScreenedNuclearRecoil::ResetTables ( )

clear precomputed screening tables

Definition at line 348 of file G4ScreenedNuclearRecoil.cc.

349 {
350 
351  std::map<G4int, G4ScreenedCoulombCrossSection*>::iterator xt=
352  crossSectionHandlers.begin();
353  for(;xt != crossSectionHandlers.end(); xt++) {
354  delete (*xt).second;
355  }
356  crossSectionHandlers.clear();
357 }
std::map< G4int, G4ScreenedCoulombCrossSection * > crossSectionHandlers

Here is the caller graph for this function:

void G4ScreenedNuclearRecoil::SetCrossSectionHardening ( G4double  fraction,
G4double  HardeningFactor 
)
inline

set the cross section boost to provide faster computation of backscattering

Parameters
fractionthe fraction of particles to have their cross section boosted.
HardeningFactorthe factor by which to boost the scattering cross section.

Definition at line 368 of file G4ScreenedNuclearRecoil.hh.

368  {
369  hardeningFraction=fraction;
370  hardeningFactor=HardeningFactor;
371  }
void G4ScreenedNuclearRecoil::SetExternalCrossSectionHandler ( G4ScreenedCoulombCrossSection cs)
inline

set a function to compute screening tables, if the user needs non-standard behavior.

Parameters
csa class which constructs the screening tables.

Definition at line 387 of file G4ScreenedNuclearRecoil.hh.

387  {
389  }
G4ScreenedCoulombCrossSection * externalCrossSectionConstructor
void G4ScreenedNuclearRecoil::SetMaxEnergyForScattering ( G4double  energy)
inline

set the upper energy beyond which this process has no cross section

This funciton is used to coordinate this process with G4MSC. Typically, G4MSC should not be allowed to operate in a range which overlaps that of this process. The criterion which is most reasonable is that the transition should be somewhere in the modestly relativistic regime (500 MeV/u for example). param energy energy per nucleon for the cutoff

Definition at line 287 of file G4ScreenedNuclearRecoil.hh.

G4double energy(const ThreeVector &p, const G4double m)
G4double processMaxEnergy
the energy per nucleon beyond which the cross section is zero, to cross over to G4MSC ...

Here is the call graph for this function:

Here is the caller graph for this function:

void G4ScreenedNuclearRecoil::SetMFPScaling ( G4double  scale)
inline

set the mean free path scaling as specified

Parameters
scalethe factor by which the default MFP will be scaled. Set to less than 1 for very thin films, typically, to sample multiple scattering, or to greater than 1 for quick simulations with a very long flight path

Definition at line 319 of file G4ScreenedNuclearRecoil.hh.

319 { MFPScale=scale; }
void G4ScreenedNuclearRecoil::SetNIELPartitionFunction ( const G4VNIELPartition part)

set the pointer to a class for paritioning energy into NIEL

part the pointer to the class.

Definition at line 371 of file G4ScreenedNuclearRecoil.cc.

373 {
376 }
const G4VNIELPartition * NIELPartitionFunction
void G4ScreenedNuclearRecoil::SetPhysicsCutoff ( G4double  energy)
inline

set the energy to which screening tables are computed. Typically, this is 10 eV or so, and not often changed.

Parameters
energythe cutoff energy

Definition at line 350 of file G4ScreenedNuclearRecoil.hh.

351  ResetTables(); }
G4double energy(const ThreeVector &p, const G4double m)
void ResetTables()
clear precomputed screening tables

Here is the call graph for this function:

void G4ScreenedNuclearRecoil::SetRecoilCutoff ( G4double  energy)
inline

set the minimum energy (per nucleon) at which recoils can be generated, and the energy (per nucleon) below which all ions are stopped.

Parameters
energyenergy per nucleon

Definition at line 341 of file G4ScreenedNuclearRecoil.hh.

341 { recoilCutoff=energy; }
G4double energy(const ThreeVector &p, const G4double m)

Here is the call graph for this function:

void G4ScreenedNuclearRecoil::SetValidCollision ( G4bool  flag)
inline

Definition at line 404 of file G4ScreenedNuclearRecoil.hh.

Here is the caller graph for this function:

Friends And Related Function Documentation

friend class G4ScreenedCollisionStage
friend

Definition at line 213 of file G4ScreenedNuclearRecoil.hh.

Member Data Documentation

G4bool G4ScreenedNuclearRecoil::avoidReactions
protected

Definition at line 430 of file G4ScreenedNuclearRecoil.hh.

std::vector<G4ScreenedCollisionStage *> G4ScreenedNuclearRecoil::collisionStages
protected

Definition at line 438 of file G4ScreenedNuclearRecoil.hh.

std::map<G4int, G4ScreenedCoulombCrossSection*> G4ScreenedNuclearRecoil::crossSectionHandlers
protected

Definition at line 440 of file G4ScreenedNuclearRecoil.hh.

G4ScreenedCoulombCrossSection* G4ScreenedNuclearRecoil::externalCrossSectionConstructor
protected

Definition at line 437 of file G4ScreenedNuclearRecoil.hh.

G4bool G4ScreenedNuclearRecoil::generateRecoils
protected

Definition at line 430 of file G4ScreenedNuclearRecoil.hh.

G4double G4ScreenedNuclearRecoil::hardeningFactor
protected

Definition at line 435 of file G4ScreenedNuclearRecoil.hh.

G4double G4ScreenedNuclearRecoil::hardeningFraction
protected

Definition at line 435 of file G4ScreenedNuclearRecoil.hh.

G4double G4ScreenedNuclearRecoil::highEnergyLimit
protected

the energy per nucleon above which the MFP is constant

Definition at line 421 of file G4ScreenedNuclearRecoil.hh.

G4double G4ScreenedNuclearRecoil::IonizingLoss
protected

Definition at line 433 of file G4ScreenedNuclearRecoil.hh.

G4CoulombKinematicsInfo G4ScreenedNuclearRecoil::kinematics
protected

Definition at line 443 of file G4ScreenedNuclearRecoil.hh.

G4double G4ScreenedNuclearRecoil::lowEnergyLimit
protected

the energy per nucleon below which the MFP is zero

Definition at line 424 of file G4ScreenedNuclearRecoil.hh.

G4double G4ScreenedNuclearRecoil::MFPScale
protected

Definition at line 434 of file G4ScreenedNuclearRecoil.hh.

G4double G4ScreenedNuclearRecoil::NIEL
protected

Definition at line 433 of file G4ScreenedNuclearRecoil.hh.

const G4VNIELPartition* G4ScreenedNuclearRecoil::NIELPartitionFunction
protected

Definition at line 444 of file G4ScreenedNuclearRecoil.hh.

G4double G4ScreenedNuclearRecoil::physicsCutoff
protected

Definition at line 431 of file G4ScreenedNuclearRecoil.hh.

G4double G4ScreenedNuclearRecoil::processMaxEnergy
protected

the energy per nucleon beyond which the cross section is zero, to cross over to G4MSC

Definition at line 428 of file G4ScreenedNuclearRecoil.hh.

G4double G4ScreenedNuclearRecoil::recoilCutoff
protected

Definition at line 431 of file G4ScreenedNuclearRecoil.hh.

G4bool G4ScreenedNuclearRecoil::registerDepositedEnergy
protected

Definition at line 432 of file G4ScreenedNuclearRecoil.hh.

G4String G4ScreenedNuclearRecoil::screeningKey
protected

Definition at line 429 of file G4ScreenedNuclearRecoil.hh.

G4bool G4ScreenedNuclearRecoil::validCollision
protected

Definition at line 442 of file G4ScreenedNuclearRecoil.hh.


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