Geant4 11.1.1
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")
 
 ~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
 
virtual G4bool IsApplicable (const G4ParticleDefinition &p) override=0
 
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, size_t idxCouple) const
 
void AddEmModel (G4int, G4VEmModel *, const G4Region *region=nullptr)
 
void SetEmModel (G4VEmModel *, G4int index=0)
 
G4int NumberOfModels () const
 
G4VEmModelEmModel (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 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 &)
 
virtual G4double GetCrossSection (const G4double, const G4MaterialCutsCouple *)
 
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
 
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 const G4VProcessGetCreatorProcess () const
 
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

void InitialiseProcess (const G4ParticleDefinition *) override
 
void StreamProcessInfo (std::ostream &outFile) const override
 
- Protected Member Functions inherited from G4VEmProcess
virtual void StreamProcessInfo (std::ostream &) const
 
virtual void InitialiseProcess (const G4ParticleDefinition *)=0
 
G4VEmModelSelectModel (G4double kinEnergy, 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)
 
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
 
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
const G4MaterialCutsCouplecurrentCouple = nullptr
 
const G4MaterialcurrentMaterial = nullptr
 
G4EmBiasingManagerbiasManager = nullptr
 
std::vector< G4double > * theEnergyOfCrossSectionMax = nullptr
 
G4double mfpKinEnergy = DBL_MAX
 
G4double preStepKinEnergy = 0.0
 
G4double preStepLogKinEnergy = LOG_EKIN_MIN
 
G4double preStepLambda = 0.0
 
G4int mainSecondaries = 1
 
G4int secID = _EM
 
G4int fluoID = _Fluorescence
 
G4int augerID = _AugerElectron
 
G4int biasID = _EM
 
G4int tripletID = _TripletElectron
 
size_t currentCoupleIndex = 0
 
size_t basedCoupleIndex = 0
 
size_t coupleIdxLambda = 0
 
size_t idxLambda = 0
 
G4bool isTheMaster = true
 
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 64 of file G4eplusAnnihilation.hh.

Constructor & Destructor Documentation

◆ G4eplusAnnihilation() [1/2]

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

Definition at line 68 of file G4eplusAnnihilation.cc.

69 : G4VEmProcess(name)
70{
71 theGamma = G4Gamma::Gamma();
72 theElectron = G4Electron::Electron();
74 SetBuildTableFlag(false);
76 SetSecondaryParticle(theGamma);
78 enableAtRestDoIt = true;
80 fEntanglementModelID = G4PhysicsModelCatalog::GetModelID("model_GammaGammaEntanglement");
81}
@ fAnnihilation
@ fEmDecreasing
static G4Electron * Electron()
Definition: G4Electron.cc:93
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
static G4int GetModelID(const G4int modelIndex)
G4int mainSecondaries
void SetBuildTableFlag(G4bool val)
void SetSecondaryParticle(const G4ParticleDefinition *p)
void SetCrossSectionType(G4CrossSectionType val)
void SetStartFromNullFlag(G4bool val)
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:363
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:410

◆ ~G4eplusAnnihilation()

G4eplusAnnihilation::~G4eplusAnnihilation ( )
overridedefault

◆ 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 123 of file G4eplusAnnihilation.cc.

126{
128
131 G4double ene(0.0);
132 G4VEmModel* model = SelectModel(ene, idx);
133
134 // define new weight for primary and secondaries
136
137 // sample secondaries
138 secParticles.clear();
139 G4double gammaCut = GetGammaEnergyCut();
141 track.GetDynamicParticle(), gammaCut);
142
143 G4int num0 = (G4int)secParticles.size();
144
145 // splitting or Russian roulette
146 if(biasManager) {
148 G4double eloss = 0.0;
150 secParticles, track, model, &fParticleChange, eloss,
151 idx, gammaCut, step.GetPostStepPoint()->GetSafety());
152 if(eloss > 0.0) {
155 }
156 }
157 }
158
159 // save secondaries
160 G4int num = (G4int)secParticles.size();
161
162 // Check that entanglement is switched on... (the following flag is
163 // set by /process/em/QuantumEntanglement).
165 // ...and that we have two gammas with both gammas' energies above
166 // gammaCut (entanglement is only programmed for e+ e- -> gamma gamma).
167 G4bool entangledgammagamma = false;
168 if (entangled) {
169 if (num == 2) {
170 entangledgammagamma = true;
171 for (const auto* p: secParticles) {
172 if (p->GetDefinition() != theGamma ||
173 p->GetKineticEnergy() < gammaCut) {
174 entangledgammagamma = false;
175 }
176 }
177 }
178 }
179
180 // Prepare a shared pointer for psossible use below. If it is used, the
181 // shared pointer is copied into the tracks through G4EntanglementAuxInfo.
182 // This ensures the clip board lasts until both tracks are destroyed.
183 std::shared_ptr<G4eplusAnnihilationEntanglementClipBoard> clipBoard;
184 if (entangledgammagamma) {
185 clipBoard = std::make_shared<G4eplusAnnihilationEntanglementClipBoard>();
186 clipBoard->SetParentParticleDefinition(track.GetDefinition());
187 }
188
189 if(num > 0) {
190
193 G4double time = track.GetGlobalTime();
194
195 for (G4int i=0; i<num; ++i) {
196 if (secParticles[i]) {
199 G4double e = dp->GetKineticEnergy();
200 G4bool good = true;
201 if(ApplyCuts()) {
202 if (p == theGamma) {
203 if (e < gammaCut) { good = false; }
204 } else if (p == theElectron) {
205 if (e < GetElectronEnergyCut()) { good = false; }
206 }
207 // added secondary if it is good
208 }
209 if (good) {
210 G4Track* t = new G4Track(dp, time, track.GetPosition());
212 if (entangledgammagamma) {
213 // entangledgammagamma is only true when there are only two gammas
214 // (See code above where entangledgammagamma is calculated.)
215 if (i == 0) { // First gamma
216 clipBoard->SetTrackA(t);
217 } else if (i == 1) { // Second gamma
218 clipBoard->SetTrackB(t);
219 }
221 (fEntanglementModelID,new G4EntanglementAuxInfo(clipBoard));
222 }
223 if (biasManager) {
224 t->SetWeight(weight * biasManager->GetWeight(i));
225 } else {
226 t->SetWeight(weight);
227 }
229
230 // define type of secondary
232 else if(i < num0) {
233 if(p == theGamma) {
235 } else {
237 }
238 } else {
240 }
241 /*
242 G4cout << "Secondary(post step) has weight " << t->GetWeight()
243 << ", Ekin= " << t->GetKineticEnergy()/MeV << " MeV "
244 << GetProcessName() << " fluoID= " << fluoID
245 << " augerID= " << augerID <<G4endl;
246 */
247 } else {
248 delete dp;
249 edep += e;
250 }
251 }
252 }
254 }
255 return &fParticleChange;
256}
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)
static G4EmParameters * Instance()
G4bool QuantumEntanglement() const
void InitializeForPostStep(const G4Track &)
void SetAuxiliaryTrackInformation(G4int id, G4VAuxiliaryTrackInformation *info) const
Definition: G4Track.cc:209
void SetWeight(G4double aValue)
const G4ThreeVector & GetPosition() const
void SetTouchableHandle(const G4TouchableHandle &apValue)
G4double GetGlobalTime() const
G4ParticleDefinition * GetDefinition() const
const G4DynamicParticle * GetDynamicParticle() const
const G4TouchableHandle & GetTouchableHandle() const
void SetCreatorModelID(const G4int id)
const G4MaterialCutsCouple * GetMaterialCutsCouple() const
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()
G4VEmModel * SelectModel(G4double kinEnergy, size_t)
G4EmBiasingManager * biasManager
G4double GetElectronEnergyCut()
std::vector< G4DynamicParticle * > secParticles
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:325

◆ AtRestGetPhysicalInteractionLength()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 96 of file G4eplusAnnihilation.cc.

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

◆ InitialiseProcess()

void G4eplusAnnihilation::InitialiseProcess ( const G4ParticleDefinition )
overrideprotectedvirtual

Implements G4VEmProcess.

Definition at line 105 of file G4eplusAnnihilation.cc.

106{
107 if(!isInitialised) {
108 isInitialised = true;
109 if(nullptr == EmModel(0)) { SetEmModel(new G4eeToTwoGammaModel()); }
112 AddEmModel(1, EmModel(0));
113 }
114}
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:746
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:753
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 89 of file G4eplusAnnihilation.cc.

90{
91 return (&p == G4Positron::Positron());
92}
static G4Positron * Positron()
Definition: G4Positron.cc:93

◆ operator=()

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

◆ ProcessDescription()

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

Reimplemented from G4VProcess.

Definition at line 260 of file G4eplusAnnihilation.cc.

261{
262 out << " Positron annihilation";
264}
void ProcessDescription(std::ostream &outFile) const override

◆ StreamProcessInfo()

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

Reimplemented from G4VEmProcess.

Definition at line 118 of file G4eplusAnnihilation.cc.

119{}

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