Geant4  10.02.p02
G4PolarizedCompton Class Reference

#include <G4PolarizedCompton.hh>

+ Inheritance diagram for G4PolarizedCompton:
+ Collaboration diagram for G4PolarizedCompton:

Public Member Functions

 G4PolarizedCompton (const G4String &processName="pol-compt", G4ProcessType type=fElectromagnetic)
 
virtual ~G4PolarizedCompton ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void PrintInfo ()
 
void SetModel (const G4String &name)
 
- Public Member Functions inherited from G4VEmProcess
 G4VEmProcess (const G4String &name, G4ProcessType type=fElectromagnetic)
 
virtual ~G4VEmProcess ()
 
virtual void ProcessDescription (std::ostream &outFile) const
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual void StartTracking (G4Track *)
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii)
 
G4double CrossSectionPerVolume (G4double kineticEnergy, const G4MaterialCutsCouple *couple)
 
G4double ComputeCrossSectionPerAtom (G4double kineticEnergy, G4double Z, G4double A=0., G4double cut=0.0)
 
G4double MeanFreePath (const G4Track &track)
 
G4double GetLambda (G4double &kinEnergy, const G4MaterialCutsCouple *couple)
 
void SetLambdaBinning (G4int nbins)
 
void SetMinKinEnergy (G4double e)
 
void SetMinKinEnergyPrim (G4double e)
 
void SetMaxKinEnergy (G4double e)
 
G4PhysicsTableLambdaTable () const
 
G4PhysicsTableLambdaTablePrim () const
 
const G4ParticleDefinitionParticle () const
 
const G4ParticleDefinitionSecondaryParticle () const
 
G4VEmModelSelectModelForMaterial (G4double kinEnergy, size_t &idxRegion) const
 
void AddEmModel (G4int, G4VEmModel *, const G4Region *region=0)
 
G4VEmModelEmModel (G4int index=1) const
 
void SetEmModel (G4VEmModel *, G4int index=1)
 
void UpdateEmModel (const G4String &, G4double, G4double)
 
G4VEmModelGetModelByIndex (G4int idx=0, G4bool ver=false) const
 
const G4ElementGetCurrentElement () const
 
const G4VEmModelGetCurrentModel () const
 
void SetCrossSectionBiasingFactor (G4double f, G4bool flag=true)
 
G4double CrossSectionBiasingFactor () const
 
void ActivateForcedInteraction (G4double length=0.0, const G4String &r="", G4bool flag=true)
 
void ActivateSecondaryBiasing (const G4String &region, G4double factor, G4double energyLimit)
 
void SetIntegral (G4bool val)
 
void SetBuildTableFlag (G4bool val)
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
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)
 
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 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

virtual void InitialiseProcess (const G4ParticleDefinition *)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
- Protected Member Functions inherited from G4VEmProcess
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
 
G4VEmModelSelectModel (G4double &kinEnergy, size_t index)
 
G4PhysicsVectorLambdaPhysicsVector (const G4MaterialCutsCouple *)
 
G4int LambdaBinning () const
 
G4double MinKinEnergy () const
 
G4double MaxKinEnergy () const
 
G4double PolarAngleLimit () const
 
G4bool IsIntegral () const
 
G4double RecalculateLambda (G4double kinEnergy, const G4MaterialCutsCouple *couple)
 
G4ParticleChangeForGammaGetParticleChange ()
 
void SetParticle (const G4ParticleDefinition *p)
 
void SetSecondaryParticle (const G4ParticleDefinition *p)
 
size_t CurrentMaterialCutsCoupleIndex () const
 
G4double GetGammaEnergyCut ()
 
G4double GetElectronEnergyCut ()
 
void SetStartFromNullFlag (G4bool val)
 
void SetSplineFlag (G4bool val)
 
const G4ElementGetTargetElement () const
 
const G4IsotopeGetTargetIsotope () const
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Private Member Functions

void CleanTable ()
 
void BuildAsymmetryTable (const G4ParticleDefinition &part)
 
G4double ComputeAsymmetry (G4double energy, const G4MaterialCutsCouple *couple, const G4ParticleDefinition &particle, G4double cut, G4double &tAsymmetry)
 
G4double ComputeSaturationFactor (const G4Track &aTrack)
 
G4PolarizedComptonoperator= (const G4PolarizedCompton &right)
 
 G4PolarizedCompton (const G4PolarizedCompton &)
 

Private Attributes

G4bool buildAsymmetryTable
 
G4bool useAsymmetryTable
 
G4bool isInitialised
 
G4int mType
 
G4PolarizedComptonModelemModel
 
G4ThreeVector targetPolarization
 

Static Private Attributes

static G4PhysicsTabletheAsymmetryTable = nullptr
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Attributes inherited from G4VEmProcess
G4ParticleChangeForGamma fParticleChange
 
- 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
 

Detailed Description

Definition at line 72 of file G4PolarizedCompton.hh.

Constructor & Destructor Documentation

G4PolarizedCompton::G4PolarizedCompton ( const G4String processName = "pol-compt",
G4ProcessType  type = fElectromagnetic 
)

Definition at line 73 of file G4PolarizedCompton.cc.

References G4Electron::Electron(), emModel, fComptonScattering, MeV, G4VEmProcess::SetBuildTableFlag(), G4VEmProcess::SetMinKinEnergyPrim(), G4VProcess::SetProcessSubType(), G4VEmProcess::SetSecondaryParticle(), G4VEmProcess::SetSplineFlag(), and G4VEmProcess::SetStartFromNullFlag().

+ Here is the call graph for this function:

G4PolarizedCompton::~G4PolarizedCompton ( )
virtual

Definition at line 93 of file G4PolarizedCompton.cc.

References CleanTable().

+ Here is the call graph for this function:

G4PolarizedCompton::G4PolarizedCompton ( const G4PolarizedCompton )
private

Member Function Documentation

void G4PolarizedCompton::BuildAsymmetryTable ( const G4ParticleDefinition part)
private

Definition at line 285 of file G4PolarizedCompton.cc.

References CleanTable(), ComputeAsymmetry(), emax, G4INCL::KinematicsUtils::energy(), G4PhysicsVector::Energy(), G4PhysicsTable::GetFlag(), G4ProductionCutsTable::GetMaterialCutsCouple(), G4ProductionCutsTable::GetProductionCutsTable(), G4ProductionCutsTable::GetTableSize(), G4VEmProcess::LambdaBinning(), G4VEmProcess::MaxKinEnergy(), G4VEmProcess::MinKinEnergy(), G4PhysicsTableHelper::PreparePhysicsTable(), G4PhysicsVector::PutValue(), G4PhysicsTableHelper::SetPhysicsVector(), G4PhysicsVector::SetSpline(), and theAsymmetryTable.

Referenced by BuildPhysicsTable().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

void G4PolarizedCompton::BuildPhysicsTable ( const G4ParticleDefinition part)
protectedvirtual

Reimplemented from G4VEmProcess.

Definition at line 270 of file G4PolarizedCompton.cc.

References BuildAsymmetryTable(), buildAsymmetryTable, G4VEmProcess::BuildPhysicsTable(), emModel, and G4VProcess::GetMasterProcess().

+ Here is the call graph for this function:

void G4PolarizedCompton::CleanTable ( )
private

Definition at line 100 of file G4PolarizedCompton.cc.

References G4PhysicsTable::clearAndDestroy(), and theAsymmetryTable.

Referenced by BuildAsymmetryTable(), and ~G4PolarizedCompton().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

G4double G4PolarizedCompton::ComputeAsymmetry ( G4double  energy,
const G4MaterialCutsCouple couple,
const G4ParticleDefinition particle,
G4double  cut,
G4double tAsymmetry 
)
private

Definition at line 330 of file G4PolarizedCompton.cc.

References G4VEmModel::CrossSection(), emModel, G4PolarizedComptonModel::SetBeamPolarization(), and G4PolarizedComptonModel::SetTargetPolarization().

Referenced by BuildAsymmetryTable().

+ Here is the call graph for this function:

+ Here is the caller graph for this function:

G4double G4PolarizedCompton::ComputeSaturationFactor ( const G4Track aTrack)
private
G4double G4PolarizedCompton::GetMeanFreePath ( const G4Track aTrack,
G4double  previousStepSize,
G4ForceCondition condition 
)
protectedvirtual

Reimplemented from G4VEmProcess.

Definition at line 155 of file G4PolarizedCompton.cc.

References ComputeSaturationFactor(), DBL_MAX, G4cout, G4endl, G4VEmProcess::GetMeanFreePath(), mm, theAsymmetryTable, useAsymmetryTable, and G4VProcess::verboseLevel.

+ Here is the call graph for this function:

void G4PolarizedCompton::InitialiseProcess ( const G4ParticleDefinition )
protectedvirtual

Implements G4VEmProcess.

Definition at line 118 of file G4PolarizedCompton.cc.

References G4VEmProcess::AddEmModel(), emModel, G4VEmProcess::EmModel(), G4EmParameters::Instance(), isInitialised, G4EmParameters::MaxKinEnergy(), G4EmParameters::MinKinEnergy(), mType, G4VEmProcess::SetEmModel(), G4VEmModel::SetHighEnergyLimit(), and G4VEmModel::SetLowEnergyLimit().

+ Here is the call graph for this function:

G4bool G4PolarizedCompton::IsApplicable ( const G4ParticleDefinition p)
virtual

Implements G4VEmProcess.

Definition at line 111 of file G4PolarizedCompton.cc.

References G4Gamma::Gamma().

+ Here is the call graph for this function:

G4PolarizedCompton& G4PolarizedCompton::operator= ( const G4PolarizedCompton right)
private
G4double G4PolarizedCompton::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
protectedvirtual

Reimplemented from G4VEmProcess.

Definition at line 173 of file G4PolarizedCompton.cc.

References ComputeSaturationFactor(), G4VProcess::currentInteractionLength, DBL_MAX, G4cout, G4endl, G4INCL::Math::max(), mm, G4VEmProcess::PostStepGetPhysicalInteractionLength(), theAsymmetryTable, G4VProcess::theNumberOfInteractionLengthLeft, useAsymmetryTable, G4VProcess::verboseLevel, and x.

+ Here is the call graph for this function:

void G4PolarizedCompton::PrintInfo ( )
virtual

Implements G4VEmProcess.

Definition at line 137 of file G4PolarizedCompton.cc.

References G4VEmProcess::EmModel(), G4cout, G4endl, and G4VEmModel::GetName().

+ Here is the call graph for this function:

void G4PolarizedCompton::SetModel ( const G4String name)

Definition at line 147 of file G4PolarizedCompton.cc.

References mType.

Member Data Documentation

G4bool G4PolarizedCompton::buildAsymmetryTable
private

Definition at line 123 of file G4PolarizedCompton.hh.

Referenced by BuildPhysicsTable().

G4PolarizedComptonModel* G4PolarizedCompton::emModel
private
G4bool G4PolarizedCompton::isInitialised
private

Definition at line 126 of file G4PolarizedCompton.hh.

Referenced by InitialiseProcess().

G4int G4PolarizedCompton::mType
private

Definition at line 127 of file G4PolarizedCompton.hh.

Referenced by InitialiseProcess(), and SetModel().

G4ThreeVector G4PolarizedCompton::targetPolarization
private

Definition at line 132 of file G4PolarizedCompton.hh.

G4PhysicsTable * G4PolarizedCompton::theAsymmetryTable = nullptr
staticprivate
G4bool G4PolarizedCompton::useAsymmetryTable
private

Definition at line 124 of file G4PolarizedCompton.hh.

Referenced by GetMeanFreePath(), and PostStepGetPhysicalInteractionLength().


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