Geant4 11.3.0
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4eplusAnnihilation Class Reference

#include <G4eplusAnnihilation.hh>

+ Inheritance diagram for G4eplusAnnihilation:

Public Member Functions

 G4eplusAnnihilation (const G4String &name="annihil")
 
 ~G4eplusAnnihilation () override
 
G4bool IsApplicable (const G4ParticleDefinition &p) final
 
G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData) override
 
G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition) override
 
void ProcessDescription (std::ostream &) const override
 
G4eplusAnnihilationoperator= (const G4eplusAnnihilation &right)=delete
 
 G4eplusAnnihilation (const G4eplusAnnihilation &)=delete
 
- Public Member Functions inherited from G4VEmProcess
 G4VEmProcess (const G4String &name, G4ProcessType type=fElectromagnetic)
 
 ~G4VEmProcess () override
 
void ProcessDescription (std::ostream &outFile) const override
 
void PreparePhysicsTable (const G4ParticleDefinition &) override
 
void BuildPhysicsTable (const G4ParticleDefinition &) override
 
void StartTracking (G4Track *) override
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &) override
 
G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false) override
 
G4bool RetrievePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii) override
 
virtual G4VEmProcessGetEmProcess (const G4String &name)
 
G4double GetLambda (G4double kinEnergy, const G4MaterialCutsCouple *couple, G4double logKinEnergy)
 
G4double GetCrossSection (const G4double kinEnergy, const G4MaterialCutsCouple *couple) override
 
G4double ComputeCrossSectionPerAtom (G4double kineticEnergy, G4double Z, G4double A=0., G4double cut=0.0)
 
G4double MeanFreePath (const G4Track &track)
 
void SetLambdaBinning (G4int nbins)
 
void SetMinKinEnergy (G4double e)
 
void SetMinKinEnergyPrim (G4double e)
 
void SetMaxKinEnergy (G4double e)
 
G4PhysicsTableLambdaTable () const
 
G4PhysicsTableLambdaTablePrim () const
 
void SetLambdaTable (G4PhysicsTable *)
 
void SetLambdaTablePrim (G4PhysicsTable *)
 
std::vector< G4double > * EnergyOfCrossSectionMax () const
 
void SetEnergyOfCrossSectionMax (std::vector< G4double > *)
 
G4CrossSectionType CrossSectionType () const
 
void SetCrossSectionType (G4CrossSectionType val)
 
const G4ParticleDefinitionParticle () const
 
const G4ParticleDefinitionSecondaryParticle () const
 
G4VEmModelSelectModelForMaterial (G4double kinEnergy, std::size_t idxCouple) const
 
void AddEmModel (G4int, G4VEmModel *, const G4Region *region=nullptr)
 
void SetEmModel (G4VEmModel *, G4int index=0)
 
G4int NumberOfModels () const
 
G4VEmModelEmModel (std::size_t index=0) const
 
const G4VEmModelGetCurrentModel () const
 
G4VEmModelGetModelByIndex (G4int idx=0, G4bool ver=false) 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 SetBuildTableFlag (G4bool val)
 
void CurrentSetup (const G4MaterialCutsCouple *, G4double energy)
 
G4bool UseBaseMaterial () const
 
void BuildLambdaTable ()
 
void StreamInfo (std::ostream &outFile, const G4ParticleDefinition &, G4bool rst=false) const
 
 G4VEmProcess (G4VEmProcess &)=delete
 
G4VEmProcessoperator= (const G4VEmProcess &right)=delete
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &aName, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
G4VDiscreteProcessoperator= (const G4VDiscreteProcess &)=delete
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4VParticleChangeAlongStepDoIt (const G4Track &, const G4Step &)
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
 
- 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
 
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 const G4VProcessGetCreatorProcess () const
 
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

void InitialiseProcess (const G4ParticleDefinition *) override
 
void StreamProcessInfo (std::ostream &outFile) const override
 
- Protected Member Functions inherited from G4VEmProcess
G4VEmModelSelectModel (G4double kinEnergy, std::size_t)
 
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
 
G4ParticleChangeForGammaGetParticleChange ()
 
void SetParticle (const G4ParticleDefinition *p)
 
void SetSecondaryParticle (const G4ParticleDefinition *p)
 
std::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
 
G4int DensityIndex (G4int idx) const
 
G4double DensityFactor (G4int idx) const
 
- 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
const G4MaterialCutsCouplecurrentCouple = nullptr
 
const G4MaterialcurrentMaterial = nullptr
 
G4EmBiasingManagerbiasManager = nullptr
 
std::vector< G4double > * theEnergyOfCrossSectionMax = nullptr
 
G4double mfpKinEnergy = DBL_MAX
 
G4double preStepKinEnergy = 0.0
 
G4double preStepLambda = 0.0
 
G4int mainSecondaries = 1
 
G4int secID = _EM
 
G4int fluoID = _Fluorescence
 
G4int augerID = _AugerElectron
 
G4int biasID = _EM
 
G4int tripletID = _TripletElectron
 
std::size_t currentCoupleIndex = 0
 
std::size_t basedCoupleIndex = 0
 
std::size_t coupleIdxLambda = 0
 
std::size_t idxLambda = 0
 
G4bool isTheMaster = false
 
G4bool baseMat = false
 
std::vector< G4DynamicParticle * > secParticles
 
G4ParticleChangeForGamma fParticleChange
 
- 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 68 of file G4eplusAnnihilation.hh.

Constructor & Destructor Documentation

◆ G4eplusAnnihilation() [1/2]

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

Definition at line 71 of file G4eplusAnnihilation.cc.

72 : G4VEmProcess(name)
73{
75 SetBuildTableFlag(false);
79 enableAtRestDoIt = true;
81 fEntanglementModelID = G4PhysicsModelCatalog::GetModelID("model_GammaGammaEntanglement");
82}
@ fAnnihilation
@ fEmDecreasing
static G4Gamma * Gamma()
Definition G4Gamma.cc:81
static G4int GetModelID(const G4int modelIndex)
G4int mainSecondaries
G4VEmProcess(const G4String &name, G4ProcessType type=fElectromagnetic)
void SetBuildTableFlag(G4bool val)
void SetSecondaryParticle(const G4ParticleDefinition *p)
void SetCrossSectionType(G4CrossSectionType val)
void SetStartFromNullFlag(G4bool val)
G4bool enableAtRestDoIt
void SetProcessSubType(G4int)

Referenced by G4eplusAnnihilation(), G4PolarizedAnnihilation::G4PolarizedAnnihilation(), and operator=().

◆ ~G4eplusAnnihilation()

G4eplusAnnihilation::~G4eplusAnnihilation ( )
override

Definition at line 86 of file G4eplusAnnihilation.cc.

87{
88 delete f2GammaAtRestModel;
89 delete f3GammaAtRestModel;
90}

◆ G4eplusAnnihilation() [2/2]

G4eplusAnnihilation::G4eplusAnnihilation ( const G4eplusAnnihilation & )
delete

Member Function Documentation

◆ AtRestDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 149 of file G4eplusAnnihilation.cc.

151{
152 // positron at rest should be killed
153 fParticleChange.InitializeForPostStep(track);
154 fParticleChange.SetProposedKineticEnergy(0.);
155 fParticleChange.ProposeTrackStatus(fStopAndKill);
156
157 auto couple = step.GetPreStepPoint()->GetMaterialCutsCouple();
158 DefineMaterial(couple);
159
160 G4double gammaCut = GetGammaEnergyCut();
161
162 // apply cuts
163 if (fApplyCuts && gammaCut > CLHEP::electron_mass_c2) {
164 fParticleChange.ProposeLocalEnergyDeposit(2*CLHEP::electron_mass_c2);
165 return &fParticleChange;
166 }
167
168 // sample secondaries
169 secParticles.clear();
170 G4double edep = 0.0;
171 if (nullptr != f3GammaAtRestModel &&
172 G4UniformRand() < currentMaterial->GetIonisation()->GetOrtoPositroniumFraction()) {
173 f3GammaAtRestModel->SampleSecondaries(secParticles, edep, couple->GetMaterial());
174 } else {
175 f2GammaAtRestModel->SampleSecondaries(secParticles, edep, couple->GetMaterial());
176 }
177
178 // define new weight for primary and secondaries
179 G4double weight = fParticleChange.GetParentWeight();
180 std::size_t num0 = secParticles.size();
181
182 // Russian roulette
183 if (nullptr != biasManager) {
184 G4int idx = couple->GetIndex();
185 if (biasManager->SecondaryBiasingRegion(idx) &&
186 !biasManager->GetDirectionalSplitting()) {
187 G4VEmModel* mod = EmModel(0);
188 G4double eloss = 0.0;
189 weight *= biasManager->ApplySecondaryBiasing(secParticles, track, mod,
190 &fParticleChange, eloss,
191 idx, gammaCut);
192 edep += eloss;
193 }
194 }
195
196 // save secondaries
197 std::size_t num = secParticles.size();
198
199 // Prepare a shared pointer only for two first gamma. If it is used, the
200 // shared pointer is copied into the tracks through G4EntanglementAuxInfo.
201 // This ensures the clip board lasts until both tracks are destroyed.
202 // It is assumed that 2 first secondary particles are the most energetic gamma
203 std::shared_ptr<G4eplusAnnihilationEntanglementClipBoard> clipBoard;
204 if (fEntangled && num >= 2) {
205 clipBoard = std::make_shared<G4eplusAnnihilationEntanglementClipBoard>();
206 clipBoard->SetParentParticleDefinition(track.GetDefinition());
207 }
208
209 if (num > 0) {
210 const G4double time = track.GetGlobalTime();
211 const G4ThreeVector& pos = track.GetPosition();
212 auto touch = track.GetTouchableHandle();
213 for (std::size_t i=0; i<num; ++i) {
214 G4DynamicParticle* dp = secParticles[i];
215 G4Track* t = new G4Track(dp, time, pos);
216 t->SetTouchableHandle(touch);
217 if (fEntangled && i < 2) {
218 // entangledgammagamma is only true when there are only two gammas
219 // (See code above where entangledgammagamma is calculated.)
220 if (nullptr != clipBoard) {
221 if (i == 0) { // First gamma
222 clipBoard->SetTrackA(t);
223 } else if (i == 1) { // Second gamma
224 clipBoard->SetTrackB(t);
225 }
226 }
228 (fEntanglementModelID, new G4EntanglementAuxInfo(clipBoard));
229 }
230 if (nullptr != biasManager) {
231 t->SetWeight(weight * biasManager->GetWeight((G4int)i));
232 } else {
233 t->SetWeight(weight);
234 }
235 pParticleChange->AddSecondary(t);
236
237 // define type of secondary
238 if (i < num0) {
240 }
241 else {
243 }
244 }
245 }
246 fParticleChange.ProposeLocalEnergyDeposit(edep);
247 return &fParticleChange;
248}
CLHEP::Hep3Vector G4ThreeVector
@ fStopAndKill
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4UniformRand()
Definition Randomize.hh:52
void SetAuxiliaryTrackInformation(G4int id, G4VAuxiliaryTrackInformation *info) const
Definition G4Track.cc:215
void SetWeight(G4double aValue)
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
G4ParticleDefinition * GetDefinition() const
const G4TouchableHandle & GetTouchableHandle() const
void SetCreatorModelID(const G4int id)
void DefineMaterial(const G4MaterialCutsCouple *couple)
G4double GetGammaEnergyCut()
G4VEmModel * EmModel(std::size_t index=0) const
G4EmBiasingManager * biasManager
std::vector< G4DynamicParticle * > secParticles
G4ParticleChangeForGamma fParticleChange
const G4Material * currentMaterial
G4VParticleChange * pParticleChange

◆ AtRestGetPhysicalInteractionLength()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 101 of file G4eplusAnnihilation.cc.

103{
105 return 0.0;
106}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced

◆ InitialiseProcess()

void G4eplusAnnihilation::InitialiseProcess ( const G4ParticleDefinition * )
overrideprotectedvirtual

Implements G4VEmProcess.

Definition at line 110 of file G4eplusAnnihilation.cc.

111{
112 if (!isInitialised) {
113 isInitialised = true;
114 if(nullptr == EmModel(0)) { SetEmModel(new G4eeToTwoGammaModel()); }
117 AddEmModel(1, EmModel(0));
118 }
119 auto param = G4EmParameters::Instance();
120
121 // AtRest models should be chosen only once
122 if (nullptr == f2GammaAtRestModel) {
123 auto type = param->PositronAtRestModelType();
124 if (type == fAllisonPositronium) {
125 f2GammaAtRestModel = new G4AllisonPositronAtRestModel();
126 } else if (type == fOrePowell) {
127 f2GammaAtRestModel = new G4AllisonPositronAtRestModel();
128 f3GammaAtRestModel = new G4OrePowellAtRestModel();
129 } else if (type == fOrePowellPolar) {
130 f2GammaAtRestModel = new G4AllisonPositronAtRestModel();
131 f3GammaAtRestModel = new G4PolarizedOrePowellAtRestModel();
132 } else {
133 f2GammaAtRestModel = new G4SimplePositronAtRestModel();
134 }
135 }
136 // Check that entanglement is switched on
137 // It may be set by the UI command "/process/em/QuantumEntanglement true".
138 fEntangled = param->QuantumEntanglement();
139 fApplyCuts = param->ApplyCuts();
140}
@ fAllisonPositronium
@ fOrePowell
@ fOrePowellPolar
static G4EmParameters * Instance()
void SetHighEnergyLimit(G4double)
void SetLowEnergyLimit(G4double)
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

Reimplemented from G4VProcess.

Definition at line 94 of file G4eplusAnnihilation.cc.

95{
96 return (&p == G4Positron::Positron());
97}
static G4Positron * Positron()
Definition G4Positron.cc:90

◆ operator=()

G4eplusAnnihilation & G4eplusAnnihilation::operator= ( const G4eplusAnnihilation & right)
delete

◆ ProcessDescription()

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

Reimplemented from G4VProcess.

Reimplemented in G4PolarizedAnnihilation.

Definition at line 252 of file G4eplusAnnihilation.cc.

253{
254 out << " Positron annihilation";
256}
void ProcessDescription(std::ostream &outFile) const override

◆ StreamProcessInfo()

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

Reimplemented from G4VEmProcess.

Definition at line 144 of file G4eplusAnnihilation.cc.

145{}

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