Geant4 10.7.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) final
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData) override
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition) override
 
virtual void ProcessDescription (std::ostream &) const override
 
- Public Member Functions inherited from G4VEmProcess
 G4VEmProcess (const G4String &name, G4ProcessType type=fElectromagnetic)
 
virtual ~G4VEmProcess ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &p) override=0
 
virtual void PrintInfo ()
 
virtual void ProcessDescription (std::ostream &outFile) const override
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &) override
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &) override
 
virtual void StartTracking (G4Track *) override
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
virtual G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &) override
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
 
virtual G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
 
virtual G4VEmProcessGetEmProcess (const G4String &name)
 
G4double CrossSectionPerVolume (G4double kineticEnergy, const G4MaterialCutsCouple *couple, G4double logKinEnergy=DBL_MAX)
 
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)
 
G4double GetLambda (G4double kinEnergy, const G4MaterialCutsCouple *couple, G4double logKinEnergy)
 
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=nullptr)
 
void SetEmModel (G4VEmModel *, G4int index=0)
 
G4VEmModelEmModel (size_t index=0) const
 
void UpdateEmModel (const G4String &, G4double, G4double)
 
G4int GetNumberOfModels () const
 
G4int GetNumberOfRegionModels (size_t couple_index) const
 
G4VEmModelGetRegionModel (G4int idx, size_t couple_index) const
 
G4VEmModelGetModelByIndex (G4int idx=0, G4bool ver=false) const
 
const G4VEmModelGetCurrentModel () const
 
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 SetEmMasterProcess (const G4VEmProcess *)
 
void SetIntegral (G4bool val)
 
void SetBuildTableFlag (G4bool val)
 
void CurrentSetup (const G4MaterialCutsCouple *, G4double energy)
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &aName, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
G4VDiscreteProcessoperator= (const G4VDiscreteProcess &)=delete
 
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 ()
 
G4VProcessoperator= (const G4VProcess &)=delete
 
G4bool operator== (const G4VProcess &right) const
 
G4bool 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
 
virtual void ProcessDescription (std::ostream &outfile) 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 StreamProcessInfo (std::ostream &outFile) const override
 
virtual void InitialiseProcess (const G4ParticleDefinition *) override
 
- Protected Member Functions inherited from G4VEmProcess
virtual void StreamProcessInfo (std::ostream &) const
 
virtual void InitialiseProcess (const G4ParticleDefinition *)=0
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
 
G4VEmModelSelectModel (G4double kinEnergy, size_t index)
 
virtual G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
G4PhysicsVectorLambdaPhysicsVector (const G4MaterialCutsCouple *)
 
void DefineMaterial (const G4MaterialCutsCouple *couple)
 
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
 
const G4MaterialCutsCoupleMaterialCutsCouple () const
 
G4bool ApplyCuts () const
 
G4double GetGammaEnergyCut ()
 
G4double GetElectronEnergyCut ()
 
void SetStartFromNullFlag (G4bool val)
 
void SetSplineFlag (G4bool val)
 
const G4ElementGetTargetElement () const
 
const G4IsotopeGetTargetIsotope () const
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double prevStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Attributes inherited from G4VEmProcess
G4LossTableManagerlManager
 
G4EmParameterstheParameters
 
G4EmBiasingManagerbiasManager
 
const G4ParticleDefinitiontheGamma
 
const G4ParticleDefinitiontheElectron
 
G4ParticleChangeForGamma fParticleChange
 
std::vector< G4DynamicParticle * > secParticles
 
const G4MaterialCutsCouplecurrentCouple
 
const G4MaterialcurrentMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t currentCoupleIndex
 
size_t basedCoupleIndex
 
G4int mainSecondaries
 
G4int secID
 
G4int fluoID
 
G4int augerID
 
G4int biasID
 
G4bool isTheMaster
 
G4double mfpKinEnergy
 
G4double preStepKinEnergy
 
G4double preStepLogKinEnergy
 
G4double preStepLambda
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4ParticleChange aParticleChange
 
G4double theNumberOfInteractionLengthLeft = -1.0
 
G4double currentInteractionLength = -1.0
 
G4double theInitialNumberOfInteractionLength = -1.0
 
G4String theProcessName
 
G4String thePhysicsTableFileName
 
G4ProcessType theProcessType = fNotDefined
 
G4int theProcessSubType = -1
 
G4double thePILfactor = 1.0
 
G4int verboseLevel = 0
 
G4bool enableAtRestDoIt = true
 
G4bool enableAlongStepDoIt = true
 
G4bool enablePostStepDoIt = true
 

Detailed Description

Definition at line 64 of file G4eplusAnnihilation.hh.

Constructor & Destructor Documentation

◆ G4eplusAnnihilation()

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

Definition at line 58 of file G4eplusAnnihilation.cc.

59 : G4VEmProcess(name), isInitialised(false)
60{
61 theGamma = G4Gamma::Gamma();
62 SetIntegral(true);
63 SetBuildTableFlag(false);
65 SetSecondaryParticle(theGamma);
67 enableAtRestDoIt = true;
69}
@ fAnnihilation
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
void SetIntegral(G4bool val)
G4int mainSecondaries
void SetBuildTableFlag(G4bool val)
void SetSecondaryParticle(const G4ParticleDefinition *p)
void SetStartFromNullFlag(G4bool val)
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:359
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:406

◆ ~G4eplusAnnihilation()

G4eplusAnnihilation::~G4eplusAnnihilation ( )
virtual

Definition at line 73 of file G4eplusAnnihilation.cc.

74{}

Member Function Documentation

◆ AtRestDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 112 of file G4eplusAnnihilation.cc.

115{
117
119 size_t idx = CurrentMaterialCutsCoupleIndex();
120 G4double ene(0.0);
121 G4VEmModel* model = SelectModel(ene, idx);
122
123 // define new weight for primary and secondaries
125
126 // sample secondaries
127 secParticles.clear();
128 G4double gammaCut = GetGammaEnergyCut();
130 track.GetDynamicParticle(), gammaCut);
131
132 G4int num0 = secParticles.size();
133
134 // splitting or Russian roulette
135 if(biasManager) {
137 G4double eloss = 0.0;
139 secParticles, track, model, &fParticleChange, eloss,
140 idx, gammaCut, step.GetPostStepPoint()->GetSafety());
141 if(eloss > 0.0) {
144 }
145 }
146 }
147 // save secondaries
148 G4int num = secParticles.size();
149 if(num > 0) {
150
153 G4double time = track.GetGlobalTime();
154
155 for (G4int i=0; i<num; ++i) {
156 if (secParticles[i]) {
159 G4double e = dp->GetKineticEnergy();
160 G4bool good = true;
161 if(ApplyCuts()) {
162 if (p == theGamma) {
163 if (e < gammaCut) { good = false; }
164 } else if (p == theElectron) {
165 if (e < GetElectronEnergyCut()) { good = false; }
166 }
167 // added secondary if it is good
168 }
169 if (good) {
170 G4Track* t = new G4Track(dp, time, track.GetPosition());
172 if (biasManager) {
173 t->SetWeight(weight * biasManager->GetWeight(i));
174 } else {
175 t->SetWeight(weight);
176 }
178
179 // define type of secondary
181 else if(i < num0) {
182 if(p == theGamma) {
184 } else {
186 }
187 } else {
189 }
190 /*
191 G4cout << "Secondary(post step) has weight " << t->GetWeight()
192 << ", Ekin= " << t->GetKineticEnergy()/MeV << " MeV "
193 << GetProcessName() << " fluoID= " << fluoID
194 << " augerID= " << augerID <<G4endl;
195 */
196 } else {
197 delete dp;
198 edep += e;
199 }
200 }
201 }
203 }
204 return &fParticleChange;
205}
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetKineticEnergy() const
G4double ApplySecondaryBiasing(std::vector< G4DynamicParticle * > &, const G4Track &track, G4VEmModel *currentModel, G4ParticleChangeForGamma *pParticleChange, G4double &eloss, G4int coupleIdx, G4double tcut, G4double safety=0.0)
G4double GetWeight(G4int i)
G4bool SecondaryBiasingRegion(G4int coupleIdx)
void InitializeForPostStep(const G4Track &)
void SetWeight(G4double aValue)
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
const G4DynamicParticle * GetDynamicParticle() const
const G4TouchableHandle & GetTouchableHandle() const
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
void SetCreatorModelIndex(G4int idx)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
void DefineMaterial(const G4MaterialCutsCouple *couple)
G4bool ApplyCuts() const
G4double GetGammaEnergyCut()
G4EmBiasingManager * biasManager
G4double GetElectronEnergyCut()
G4VEmModel * SelectModel(G4double kinEnergy, size_t index)
std::vector< G4DynamicParticle * > secParticles
const G4ParticleDefinition * theElectron
const G4MaterialCutsCouple * MaterialCutsCouple() const
G4ParticleChangeForGamma fParticleChange
size_t CurrentMaterialCutsCoupleIndex() const
G4double GetParentWeight() const
G4double GetLocalEnergyDeposit() const
void AddSecondary(G4Track *aSecondary)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
void SetNumberOfSecondaries(G4int totSecondaries)
G4VParticleChange * pParticleChange
Definition: G4VProcess.hh:321

◆ AtRestGetPhysicalInteractionLength()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 85 of file G4eplusAnnihilation.cc.

87{
89 return 0.0;
90}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced

◆ InitialiseProcess()

void G4eplusAnnihilation::InitialiseProcess ( const G4ParticleDefinition )
overrideprotectedvirtual

Implements G4VEmProcess.

Definition at line 94 of file G4eplusAnnihilation.cc.

95{
96 if(!isInitialised) {
97 isInitialised = true;
98 if(!EmModel(0)) { SetEmModel(new G4eeToTwoGammaModel()); }
101 AddEmModel(1, EmModel(0));
102 }
103}
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:757
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:764
G4VEmModel * EmModel(size_t index=0) const
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=nullptr)
void SetEmModel(G4VEmModel *, G4int index=0)
G4double MaxKinEnergy() const
G4double MinKinEnergy() const

◆ IsApplicable()

G4bool G4eplusAnnihilation::IsApplicable ( const G4ParticleDefinition p)
finalvirtual

Implements G4VEmProcess.

Definition at line 78 of file G4eplusAnnihilation.cc.

79{
80 return (&p == G4Positron::Positron());
81}
static G4Positron * Positron()
Definition: G4Positron.cc:93

◆ ProcessDescription()

void G4eplusAnnihilation::ProcessDescription ( std::ostream &  out) const
overridevirtual

Reimplemented from G4VEmProcess.

Definition at line 209 of file G4eplusAnnihilation.cc.

210{
211 out << " Positron annihilation";
213}
virtual void ProcessDescription(std::ostream &outFile) const override

◆ StreamProcessInfo()

void G4eplusAnnihilation::StreamProcessInfo ( std::ostream &  outFile) const
overrideprotectedvirtual

Reimplemented from G4VEmProcess.

Definition at line 107 of file G4eplusAnnihilation.cc.

108{}

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