Geant4  10.00.p01
G4PenelopeCrossSection Class Reference

#include <G4PenelopeCrossSection.hh>

+ Collaboration diagram for G4PenelopeCrossSection:

Public Member Functions

 G4PenelopeCrossSection (size_t nOfEnergyPoints, size_t nOfShells=0)
 
 ~G4PenelopeCrossSection ()
 
G4double GetTotalCrossSection (G4double energy) const
 Returns total cross section at the given energy. More...
 
G4double GetHardCrossSection (G4double energy) const
 Returns hard cross section at the given energy. More...
 
G4double GetSoftStoppingPower (G4double energy) const
 Returns the total stopping power due to soft collisions. More...
 
G4double GetShellCrossSection (size_t shellID, G4double energy) const
 Returns the hard cross section for the given shell (per molecule) More...
 
G4double GetNormalizedShellCrossSection (size_t shellID, G4double energy) const
 Returns the hard cross section for the given shell (normalized to 1) More...
 
size_t GetNumberOfShells () const
 
void AddCrossSectionPoint (size_t binNumber, G4double energy, G4double XH0, G4double XH1, G4double XH2, G4double XS0, G4double XS1, G4double XS2)
 Public interface for the master thread. More...
 
void AddShellCrossSectionPoint (size_t binNumber, size_t shellID, G4double energy, G4double xs)
 
void NormalizeShellCrossSections ()
 

Private Member Functions

G4PenelopeCrossSectionoperator= (const G4PenelopeCrossSection &right)
 
 G4PenelopeCrossSection (const G4PenelopeCrossSection &)
 

Private Attributes

G4bool isNormalized
 
size_t numberOfEnergyPoints
 
size_t numberOfShells
 
G4PhysicsTablesoftCrossSections
 
G4PhysicsTablehardCrossSections
 
G4PhysicsTableshellCrossSections
 
G4PhysicsTableshellNormalizedCrossSections
 

Detailed Description

Definition at line 72 of file G4PenelopeCrossSection.hh.

Constructor & Destructor Documentation

G4PenelopeCrossSection::G4PenelopeCrossSection ( size_t  nOfEnergyPoints,
size_t  nOfShells = 0 
)

Definition at line 44 of file G4PenelopeCrossSection.cc.

References FatalException, G4endl, G4Exception(), hardCrossSections, isNormalized, numberOfEnergyPoints, numberOfShells, G4PhysicsTable::push_back(), shellCrossSections, shellNormalizedCrossSections, and softCrossSections.

+ Here is the call graph for this function:

G4PenelopeCrossSection::~G4PenelopeCrossSection ( )
G4PenelopeCrossSection::G4PenelopeCrossSection ( const G4PenelopeCrossSection )
private

Member Function Documentation

void G4PenelopeCrossSection::AddCrossSectionPoint ( size_t  binNumber,
G4double  energy,
G4double  XH0,
G4double  XH1,
G4double  XH2,
G4double  XS0,
G4double  XS1,
G4double  XS2 
)

Public interface for the master thread.

Definition at line 123 of file G4PenelopeCrossSection.cc.

References cm2, eV, G4cout, G4endl, hardCrossSections, G4INCL::Math::max(), numberOfEnergyPoints, G4PhysicsFreeVector::PutValue(), and softCrossSections.

Referenced by G4PenelopeIonisationXSHandler::BuildXSTable(), and G4PenelopeBremsstrahlungModel::BuildXSTable().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void G4PenelopeCrossSection::AddShellCrossSectionPoint ( size_t  binNumber,
size_t  shellID,
G4double  energy,
G4double  xs 
)

Definition at line 183 of file G4PenelopeCrossSection.cc.

References cm2, G4cout, G4endl, G4INCL::Math::max(), numberOfEnergyPoints, numberOfShells, G4PhysicsFreeVector::PutValue(), and shellCrossSections.

Referenced by G4PenelopeIonisationXSHandler::BuildXSTable().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

G4double G4PenelopeCrossSection::GetHardCrossSection ( G4double  energy) const

Returns hard cross section at the given energy.

Definition at line 268 of file G4PenelopeCrossSection.cc.

References G4cout, G4endl, G4PhysicsVector::GetVectorLength(), hardCrossSections, numberOfEnergyPoints, and G4PhysicsVector::Value().

Referenced by G4PenelopeBremsstrahlungModel::CrossSectionPerVolume(), and G4PenelopeIonisationModel::CrossSectionPerVolume().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

G4double G4PenelopeCrossSection::GetNormalizedShellCrossSection ( size_t  shellID,
G4double  energy 
) const

Returns the hard cross section for the given shell (normalized to 1)

Definition at line 363 of file G4PenelopeCrossSection.cc.

References G4cout, G4endl, G4PhysicsVector::GetVectorLength(), isNormalized, numberOfEnergyPoints, numberOfShells, shellNormalizedCrossSections, and G4PhysicsVector::Value().

Referenced by G4PenelopeIonisationModel::SampleFinalStateElectron(), and G4PenelopeIonisationModel::SampleFinalStatePositron().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

size_t G4PenelopeCrossSection::GetNumberOfShells ( ) const
inline

Definition at line 94 of file G4PenelopeCrossSection.hh.

References numberOfShells.

G4double G4PenelopeCrossSection::GetShellCrossSection ( size_t  shellID,
G4double  energy 
) const

Returns the hard cross section for the given shell (per molecule)

Definition at line 327 of file G4PenelopeCrossSection.cc.

References G4cout, G4endl, G4PhysicsVector::GetVectorLength(), numberOfEnergyPoints, numberOfShells, shellCrossSections, and G4PhysicsVector::Value().

Referenced by G4PenelopeIonisationCrossSection::CrossSection().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

G4double G4PenelopeCrossSection::GetSoftStoppingPower ( G4double  energy) const

Returns the total stopping power due to soft collisions.

Definition at line 298 of file G4PenelopeCrossSection.cc.

References G4cout, G4endl, G4PhysicsVector::GetVectorLength(), numberOfEnergyPoints, softCrossSections, and G4PhysicsVector::Value().

Referenced by G4PenelopeBremsstrahlungModel::ComputeDEDXPerVolume(), and G4PenelopeIonisationModel::ComputeDEDXPerVolume().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

G4double G4PenelopeCrossSection::GetTotalCrossSection ( G4double  energy) const

Returns total cross section at the given energy.

Definition at line 224 of file G4PenelopeCrossSection.cc.

References G4cout, G4endl, G4PhysicsVector::GetVectorLength(), hardCrossSections, numberOfEnergyPoints, softCrossSections, and G4PhysicsVector::Value().

Referenced by G4PenelopeIonisationModel::CrossSectionPerVolume().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void G4PenelopeCrossSection::NormalizeShellCrossSections ( )

Definition at line 411 of file G4PenelopeCrossSection.cc.

References G4cout, G4endl, G4PhysicsVector::GetLowEdgeEnergy(), isNormalized, numberOfEnergyPoints, numberOfShells, G4PhysicsFreeVector::PutValue(), shellCrossSections, and shellNormalizedCrossSections.

Referenced by G4PenelopeIonisationXSHandler::BuildXSTable().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

G4PenelopeCrossSection& G4PenelopeCrossSection::operator= ( const G4PenelopeCrossSection right)
private

Member Data Documentation

G4PhysicsTable* G4PenelopeCrossSection::hardCrossSections
private
G4bool G4PenelopeCrossSection::isNormalized
private
G4PhysicsTable* G4PenelopeCrossSection::shellNormalizedCrossSections
private
G4PhysicsTable* G4PenelopeCrossSection::softCrossSections
private

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