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

#include <G4eplusPolarizedAnnihilation.hh>

+ Inheritance diagram for G4eplusPolarizedAnnihilation:

Public Member Functions

 G4eplusPolarizedAnnihilation (const G4String &name="pol-annihil")
 
virtual ~G4eplusPolarizedAnnihilation ()
 
virtual void PrintInfo () override
 
virtual G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &) override
 
- Public Member Functions inherited from G4eplusAnnihilation
 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 &)
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Member Functions inherited from G4eplusAnnihilation
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 ()
 
- 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 62 of file G4eplusPolarizedAnnihilation.hh.

Constructor & Destructor Documentation

◆ G4eplusPolarizedAnnihilation()

G4eplusPolarizedAnnihilation::G4eplusPolarizedAnnihilation ( const G4String name = "pol-annihil")
explicit

Definition at line 73 of file G4eplusPolarizedAnnihilation.cc.

74 : G4eplusAnnihilation(name), isInitialised(false),
75 theAsymmetryTable(nullptr),
76 theTransverseAsymmetryTable(nullptr)
77{
78 emModel = new G4PolarizedAnnihilationModel();
79 SetEmModel(emModel);
80}
void SetEmModel(G4VEmModel *, G4int index=0)

◆ ~G4eplusPolarizedAnnihilation()

G4eplusPolarizedAnnihilation::~G4eplusPolarizedAnnihilation ( )
virtual

Definition at line 84 of file G4eplusPolarizedAnnihilation.cc.

85{
86 CleanTables();
87}

Member Function Documentation

◆ BuildPhysicsTable()

void G4eplusPolarizedAnnihilation::BuildPhysicsTable ( const G4ParticleDefinition part)
overridevirtual

Reimplemented from G4VEmProcess.

Definition at line 235 of file G4eplusPolarizedAnnihilation.cc.

237{
239 G4bool isMaster = true;
240 const G4eplusPolarizedAnnihilation* masterProcess =
241 static_cast<const G4eplusPolarizedAnnihilation*>(GetMasterProcess());
242 if(masterProcess && masterProcess != this) { isMaster = false; }
243 if(isMaster) { BuildAsymmetryTables(part); }
244}
bool G4bool
Definition: G4Types.hh:86
virtual void BuildPhysicsTable(const G4ParticleDefinition &) override
const G4VProcess * GetMasterProcess() const
Definition: G4VProcess.hh:518

◆ GetMeanFreePath()

G4double G4eplusPolarizedAnnihilation::GetMeanFreePath ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
overridevirtual

Reimplemented from G4VEmProcess.

Definition at line 107 of file G4eplusPolarizedAnnihilation.cc.

110{
111 G4double mfp = G4VEmProcess::GetMeanFreePath(track, previousStepSize, condition);
112
113 if(theAsymmetryTable && theTransverseAsymmetryTable && mfp < DBL_MAX) {
114 mfp *= ComputeSaturationFactor(track);
115 }
116 if (verboseLevel>=2) {
117 G4cout << "G4eplusPolarizedAnnihilation::MeanFreePath: "
118 << mfp / mm << " mm " << G4endl;
119 }
120 return mfp;
121}
G4double condition(const G4ErrorSymMatrix &m)
double G4double
Definition: G4Types.hh:83
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
virtual G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4int verboseLevel
Definition: G4VProcess.hh:356
#define DBL_MAX
Definition: templates.hh:62

◆ PostStepGetPhysicalInteractionLength()

G4double G4eplusPolarizedAnnihilation::PostStepGetPhysicalInteractionLength ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
)
overridevirtual

Reimplemented from G4VEmProcess.

Definition at line 125 of file G4eplusPolarizedAnnihilation.cc.

129{
130 // save previous values
133
134 // *** compute unpolarized step limit ***
135 // this changes theNumberOfInteractionLengthLeft and currentInteractionLength
137 previousStepSize,
138 condition);
139 G4double x0 = x;
140 G4double satFact = 1.0;
141
142 // *** add corrections on polarisation ***
143 if(theAsymmetryTable && theTransverseAsymmetryTable && x < DBL_MAX) {
144 satFact = ComputeSaturationFactor(track);
145 G4double curLength = currentInteractionLength*satFact;
146 G4double prvLength = iLength*satFact;
147 if(nLength > 0.0) {
149 std::max(nLength - previousStepSize/prvLength, 0.0);
150 }
151 x = theNumberOfInteractionLengthLeft * curLength;
152 }
153 if (verboseLevel>=2) {
154 G4cout << "G4eplusPolarizedAnnihilation::PostStepGPIL: "
155 << std::setprecision(8) << x/mm << " mm;" << G4endl
156 << " unpolarized value: "
157 << std::setprecision(8) << x0/mm << " mm." << G4endl;
158 }
159 return x;
160}
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition) override
G4double currentInteractionLength
Definition: G4VProcess.hh:335
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:331

◆ PrintInfo()

void G4eplusPolarizedAnnihilation::PrintInfo ( )
overridevirtual

Reimplemented from G4VEmProcess.

Definition at line 329 of file G4eplusPolarizedAnnihilation.cc.

330{
331 G4cout << " Polarized model for annihilation into 2 photons"
332 << G4endl;
333}

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