Geant4  10.02.p03
G4eplusAnnihilation Class Reference

#include <G4eplusAnnihilation.hh>

Inheritance diagram for G4eplusAnnihilation:
Collaboration diagram for G4eplusAnnihilation:

Public Member Functions

 G4eplusAnnihilation (const G4String &name="annihil")
 
virtual ~G4eplusAnnihilation ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &p)
 
virtual G4VParticleChange * AtRestDoIt (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 void ProcessDescription (std::ostream &outFile) const
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
virtual void StartTracking (G4Track *)
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4VParticleChange * PostStepDoIt (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 G4VParticleChange * AlongStepDoIt (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 *)
 
- Protected Member Functions inherited from G4VEmProcess
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
 
G4VEmModelSelectModel (G4double &kinEnergy, size_t index)
 
virtual G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
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)
 
G4ParticleChangeForGamma * GetParticleChange ()
 
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 Attributes

G4bool isInitialised
 
const G4ParticleDefinitiontheGamma
 

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
 
G4VParticleChange * pParticleChange
 
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 65 of file G4eplusAnnihilation.hh.

Constructor & Destructor Documentation

◆ G4eplusAnnihilation()

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

Definition at line 67 of file G4eplusAnnihilation.cc.

68  : G4VEmProcess(name), isInitialised(false)
69 {
71  SetIntegral(true);
72  SetBuildTableFlag(false);
73  SetStartFromNullFlag(false);
76  enableAtRestDoIt = true;
77 }
G4VEmProcess(const G4String &name, G4ProcessType type=fElectromagnetic)
Definition: G4VEmProcess.cc:91
void SetBuildTableFlag(G4bool val)
void SetStartFromNullFlag(G4bool val)
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:350
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:432
void SetIntegral(G4bool val)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
const G4ParticleDefinition * theGamma
void SetSecondaryParticle(const G4ParticleDefinition *p)
Here is the call graph for this function:

◆ ~G4eplusAnnihilation()

G4eplusAnnihilation::~G4eplusAnnihilation ( )
virtual

Definition at line 81 of file G4eplusAnnihilation.cc.

82 {}

Member Function Documentation

◆ AtRestDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 120 of file G4eplusAnnihilation.cc.

129 {
130  fParticleChange.InitializeForPostStep(aTrack);
131 
132  G4double cosTeta = 2.*G4UniformRand()-1.;
133  G4double sinTeta = sqrt((1.-cosTeta)*(1.0 + cosTeta));
134  G4double phi = twopi * G4UniformRand();
135  G4ThreeVector dir(sinTeta*cos(phi), sinTeta*sin(phi), cosTeta);
136  phi = twopi * G4UniformRand();
137  G4ThreeVector pol(cos(phi), sin(phi), 0.0);
138  pol.rotateUz(dir);
139 
140  // e+ parameters
141  G4double weight = aTrack.GetWeight();
142  G4double time = aTrack.GetGlobalTime();
143 
144  // add gammas
145  fParticleChange.SetNumberOfSecondaries(2);
146  G4DynamicParticle* dp =
148  dp->SetPolarization(pol.x(),pol.y(),pol.z());
149  G4Track* track = new G4Track(dp, time, aTrack.GetPosition());
150  track->SetTouchableHandle(aTrack.GetTouchableHandle());
151  track->SetWeight(weight);
152  pParticleChange->AddSecondary(track);
153 
155  dp->SetPolarization(-pol.x(),-pol.y(),-pol.z());
156  track = new G4Track(dp, time, aTrack.GetPosition());
157  track->SetTouchableHandle(aTrack.GetTouchableHandle());
158  track->SetWeight(weight);
159  pParticleChange->AddSecondary(track);
160 
161  // Kill the incident positron
162  //
163  fParticleChange.ProposeTrackStatus(fStopAndKill);
164  return &fParticleChange;
165 }
G4ParticleChangeForGamma fParticleChange
TDirectory * dir
double weight
Definition: plottest35.C:25
#define G4UniformRand()
Definition: Randomize.hh:97
static const double twopi
Definition: G4SIunits.hh:75
float electron_mass_c2
Definition: hepunit.py:274
void SetPolarization(G4double polX, G4double polY, G4double polZ)
const G4ParticleDefinition * theGamma
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:283
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:

◆ AtRestGetPhysicalInteractionLength()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 93 of file G4eplusAnnihilation.cc.

95 {
96  *condition = NotForced;
97  return 0.0;
98 }
G4double condition(const G4ErrorSymMatrix &m)

◆ InitialiseProcess()

void G4eplusAnnihilation::InitialiseProcess ( const G4ParticleDefinition )
protectedvirtual

Implements G4VEmProcess.

Definition at line 102 of file G4eplusAnnihilation.cc.

103 {
104  if(!isInitialised) {
105  isInitialised = true;
106  if(!EmModel(1)) { SetEmModel(new G4eeToTwoGammaModel(),1); }
109  AddEmModel(1, EmModel(1));
110  }
111 }
G4double MinKinEnergy() const
G4double MaxKinEnergy() const
G4VEmModel * EmModel(G4int index=1) const
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:725
void SetEmModel(G4VEmModel *, G4int index=1)
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=0)
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:732
Here is the call graph for this function:

◆ IsApplicable()

G4bool G4eplusAnnihilation::IsApplicable ( const G4ParticleDefinition p)
virtual

Implements G4VEmProcess.

Definition at line 86 of file G4eplusAnnihilation.cc.

87 {
88  return (&p == G4Positron::Positron());
89 }
static G4Positron * Positron()
Definition: G4Positron.cc:94
Here is the call graph for this function:

◆ PrintInfo()

void G4eplusAnnihilation::PrintInfo ( )
virtual

Implements G4VEmProcess.

Reimplemented in G4eplusPolarizedAnnihilation.

Definition at line 115 of file G4eplusAnnihilation.cc.

116 {}

Member Data Documentation

◆ isInitialised

G4bool G4eplusAnnihilation::isInitialised
private

Definition at line 94 of file G4eplusAnnihilation.hh.

◆ theGamma

const G4ParticleDefinition* G4eplusAnnihilation::theGamma
private

Definition at line 95 of file G4eplusAnnihilation.hh.


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