Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4eplusAnnihilation Class Reference

#include <G4eplusAnnihilation.hh>

+ Inheritance diagram for G4eplusAnnihilation:

Public Member Functions

 G4eplusAnnihilation (const G4String &name="annihil")
 
virtual ~G4eplusAnnihilation ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &p)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)
 
virtual void PrintInfo ()
 
- Public Member Functions inherited from G4VEmProcess
 G4VEmProcess (const G4String &name, G4ProcessType type=fElectromagnetic)
 
virtual ~G4VEmProcess ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &p)=0
 
virtual void PrintInfo ()=0
 
void PreparePhysicsTable (const G4ParticleDefinition &)
 
void BuildPhysicsTable (const G4ParticleDefinition &)
 
void PrintInfoDefinition ()
 
void StartTracking (G4Track *)
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
 
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)
 
G4int LambdaBinning () const
 
void SetMinKinEnergy (G4double e)
 
G4double MinKinEnergy () const
 
void SetMaxKinEnergy (G4double e)
 
G4double MaxKinEnergy () const
 
void SetMinKinEnergyPrim (G4double e)
 
const G4PhysicsTableLambdaTable () const
 
const G4ParticleDefinitionParticle () const
 
const G4ParticleDefinitionSecondaryParticle () const
 
G4VEmModelSelectModelForMaterial (G4double kinEnergy, size_t &idxRegion) const
 
void AddEmModel (G4int, G4VEmModel *, const G4Region *region=0)
 
void SetModel (G4VEmModel *, G4int index=1)
 
G4VEmModelModel (G4int index=1)
 
G4VEmModelEmModel (G4int index=1)
 
void SetEmModel (G4VEmModel *, G4int index=1)
 
void UpdateEmModel (const G4String &, G4double, G4double)
 
G4VEmModelGetModelByIndex (G4int idx=0, G4bool ver=false)
 
const G4ElementGetCurrentElement () 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 SetPolarAngleLimit (G4double a)
 
G4double PolarAngleLimit () const
 
void SetLambdaFactor (G4double val)
 
void SetIntegral (G4bool val)
 
G4bool IsIntegral () const
 
void SetApplyCuts (G4bool val)
 
void SetBuildTableFlag (G4bool val)
 
- 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 G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
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
 
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)=0
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &proposedSafety, G4GPILSelection *selection)=0
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)=0
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
 
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 G4bool IsApplicable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
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
 

Protected Member Functions

virtual void InitialiseProcess (const G4ParticleDefinition *)
 
- Protected Member Functions inherited from G4VEmProcess
virtual void InitialiseProcess (const G4ParticleDefinition *)=0
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
 
G4VEmModelSelectModel (G4double &kinEnergy, size_t index)
 
G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
G4PhysicsVectorLambdaPhysicsVector (const G4MaterialCutsCouple *)
 
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)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

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 67 of file G4eplusAnnihilation.hh.

Constructor & Destructor Documentation

◆ G4eplusAnnihilation()

G4eplusAnnihilation::G4eplusAnnihilation ( const G4String name = "annihil")

Definition at line 68 of file G4eplusAnnihilation.cc.

69 : G4VEmProcess(name), isInitialised(false)
70{
71 theGamma = G4Gamma::Gamma();
72 SetIntegral(true);
73 SetBuildTableFlag(false);
75 SetSecondaryParticle(theGamma);
77 enableAtRestDoIt = true;
78}
@ fAnnihilation
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void SetIntegral(G4bool val)
void SetBuildTableFlag(G4bool val)
void SetSecondaryParticle(const G4ParticleDefinition *p)
void SetStartFromNullFlag(G4bool val)
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:350
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403

◆ ~G4eplusAnnihilation()

G4eplusAnnihilation::~G4eplusAnnihilation ( )
virtual

Definition at line 82 of file G4eplusAnnihilation.cc.

83{}

Member Function Documentation

◆ AtRestDoIt()

G4VParticleChange * G4eplusAnnihilation::AtRestDoIt ( const G4Track track,
const G4Step stepData 
)
virtual

Reimplemented from G4VDiscreteProcess.

Definition at line 121 of file G4eplusAnnihilation.cc.

130{
132
133 G4double cosTeta = 2.*G4UniformRand()-1.;
134 G4double sinTeta = sqrt((1.-cosTeta)*(1.0 + cosTeta));
135 G4double phi = twopi * G4UniformRand();
136 G4ThreeVector dir(sinTeta*cos(phi), sinTeta*sin(phi), cosTeta);
137 phi = twopi * G4UniformRand();
138 G4ThreeVector pol(cos(phi), sin(phi), 0.0);
139 pol.rotateUz(dir);
140
141 fParticleChange.ProposeWeight(aTrack.GetWeight());
142 // add gammas
144 G4DynamicParticle* dp =
145 new G4DynamicParticle(theGamma, dir, electron_mass_c2);
146 dp->SetPolarization(pol.x(),pol.y(),pol.z());
148
149 dp = new G4DynamicParticle(theGamma,-dir, electron_mass_c2);
150 dp->SetPolarization(-pol.x(),-pol.y(),-pol.z());
152
153 // Kill the incident positron
154 //
156 return &fParticleChange;
157}
@ fStopAndKill
double G4double
Definition: G4Types.hh:64
#define G4UniformRand()
Definition: Randomize.hh:53
void SetPolarization(G4double polX, G4double polY, G4double polZ)
void InitializeForPostStep(const G4Track &)
void AddSecondary(G4DynamicParticle *aParticle)
G4ParticleChangeForGamma fParticleChange
void ProposeTrackStatus(G4TrackStatus status)
void ProposeWeight(G4double finalWeight)
void SetNumberOfSecondaries(G4int totSecondaries)

◆ AtRestGetPhysicalInteractionLength()

G4double G4eplusAnnihilation::AtRestGetPhysicalInteractionLength ( const G4Track track,
G4ForceCondition condition 
)
virtual

Reimplemented from G4VDiscreteProcess.

Definition at line 94 of file G4eplusAnnihilation.cc.

96{
98 return 0.0;
99}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced

◆ InitialiseProcess()

void G4eplusAnnihilation::InitialiseProcess ( const G4ParticleDefinition )
protectedvirtual

Implements G4VEmProcess.

Definition at line 103 of file G4eplusAnnihilation.cc.

104{
105 if(!isInitialised) {
106 isInitialised = true;
107 if(!EmModel(1)) { SetEmModel(new G4eeToTwoGammaModel(),1); }
110 AddEmModel(1, EmModel(1));
111 }
112}
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:592
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=0)
void SetEmModel(G4VEmModel *, G4int index=1)
G4VEmModel * EmModel(G4int index=1)
G4double MaxKinEnergy() const
G4double MinKinEnergy() const

◆ IsApplicable()

G4bool G4eplusAnnihilation::IsApplicable ( const G4ParticleDefinition p)
virtual

Implements G4VEmProcess.

Definition at line 87 of file G4eplusAnnihilation.cc.

88{
89 return (&p == G4Positron::Positron());
90}
static G4Positron * Positron()
Definition: G4Positron.cc:94

◆ PrintInfo()

void G4eplusAnnihilation::PrintInfo ( )
virtual

Implements G4VEmProcess.

Definition at line 116 of file G4eplusAnnihilation.cc.

117{}

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