Geant4 9.6.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 G4bool IsApplicable (const G4ParticleDefinition &p)
 
virtual G4VParticleChangeAtRestDoIt (const G4Track &track, const G4Step &stepData)
 
G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)
 
virtual void PrintInfo ()
 
G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
void BuildAsymmetryTable (const G4ParticleDefinition &part)
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
G4double ComputeAsymmetry (G4double energy, const G4MaterialCutsCouple *couple, const G4ParticleDefinition &particle, G4double cut, G4double &tasm)
 
- Public Member Functions inherited from G4VEmProcess
 G4VEmProcess (const G4String &name, G4ProcessType type=fElectromagnetic)
 
virtual ~G4VEmProcess ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &p)=0
 
virtual void PrintInfo ()=0
 
void PreparePhysicsTable (const G4ParticleDefinition &)
 
void BuildPhysicsTable (const G4ParticleDefinition &)
 
void PrintInfoDefinition ()
 
void StartTracking (G4Track *)
 
G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
G4VParticleChangePostStepDoIt (const G4Track &, const G4Step &)
 
G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
 
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)
 
G4int LambdaBinning () const
 
void SetMinKinEnergy (G4double e)
 
G4double MinKinEnergy () const
 
void SetMaxKinEnergy (G4double e)
 
G4double MaxKinEnergy () const
 
void SetMinKinEnergyPrim (G4double e)
 
const G4PhysicsTableLambdaTable () const
 
const G4ParticleDefinitionParticle () const
 
const G4ParticleDefinitionSecondaryParticle () const
 
G4VEmModelSelectModelForMaterial (G4double kinEnergy, size_t &idxRegion) const
 
void AddEmModel (G4int, G4VEmModel *, const G4Region *region=0)
 
void SetModel (G4VEmModel *, G4int index=1)
 
G4VEmModelModel (G4int index=1)
 
G4VEmModelEmModel (G4int index=1)
 
void SetEmModel (G4VEmModel *, G4int index=1)
 
void UpdateEmModel (const G4String &, G4double, G4double)
 
G4VEmModelGetModelByIndex (G4int idx=0, G4bool ver=false)
 
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 SetPolarAngleLimit (G4double a)
 
G4double PolarAngleLimit () const
 
void SetLambdaFactor (G4double val)
 
void SetIntegral (G4bool val)
 
G4bool IsIntegral () const
 
void SetApplyCuts (G4bool val)
 
void SetBuildTableFlag (G4bool val)
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
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 ()
 
G4int operator== (const G4VProcess &right) const
 
G4int 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
 
void SetVerboseLevel (G4int value)
 
G4int GetVerboseLevel () const
 

Protected Member Functions

virtual void InitialiseProcess (const G4ParticleDefinition *)
 
- Protected Member Functions inherited from G4VEmProcess
virtual void InitialiseProcess (const G4ParticleDefinition *)=0
 
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
 
G4VEmModelSelectModel (G4double &kinEnergy, size_t index)
 
G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
G4PhysicsVectorLambdaPhysicsVector (const G4MaterialCutsCouple *)
 
G4double RecalculateLambda (G4double kinEnergy, const G4MaterialCutsCouple *couple)
 
G4ParticleChangeForGammaGetParticleChange ()
 
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)
 
virtual G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)=0
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 

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
 
G4VParticleChangepParticleChange
 
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 63 of file G4eplusPolarizedAnnihilation.hh.

Constructor & Destructor Documentation

◆ G4eplusPolarizedAnnihilation()

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

Definition at line 75 of file G4eplusPolarizedAnnihilation.cc.

76 : G4VEmProcess(name), isInitialised(false),
77 theAsymmetryTable(NULL),
78 theTransverseAsymmetryTable(NULL)
79{
80 enableAtRestDoIt = true;
82 emModel = 0;
83}
@ fAnnihilation
G4bool enableAtRestDoIt
Definition: G4VProcess.hh:350
void SetProcessSubType(G4int)
Definition: G4VProcess.hh:403

◆ ~G4eplusPolarizedAnnihilation()

G4eplusPolarizedAnnihilation::~G4eplusPolarizedAnnihilation ( )
virtual

Definition at line 87 of file G4eplusPolarizedAnnihilation.cc.

88{
89 if (theAsymmetryTable) {
90 theAsymmetryTable->clearAndDestroy();
91 delete theAsymmetryTable;
92 }
93 if (theTransverseAsymmetryTable) {
94 theTransverseAsymmetryTable->clearAndDestroy();
95 delete theTransverseAsymmetryTable;
96 }
97}
void clearAndDestroy()

Member Function Documentation

◆ AtRestDoIt()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 359 of file G4eplusPolarizedAnnihilation.cc.

368{
370
372
373 G4double cosTeta = 2.*G4UniformRand()-1. , sinTeta = std::sqrt(1.-cosTeta*cosTeta);
374 G4double phi = twopi * G4UniformRand();
375 G4ThreeVector direction (sinTeta*std::cos(phi), sinTeta*std::sin(phi), cosTeta);
377 direction, electron_mass_c2) );
379 -direction, electron_mass_c2) );
380 // Kill the incident positron
381 //
383 return &fParticleChange;
384}
@ fStopAndKill
double G4double
Definition: G4Types.hh:64
#define G4UniformRand()
Definition: Randomize.hh:53
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
void InitializeForPostStep(const G4Track &)
void AddSecondary(G4DynamicParticle *aParticle)
G4ParticleChangeForGamma fParticleChange
void ProposeTrackStatus(G4TrackStatus status)
void SetNumberOfSecondaries(G4int totSecondaries)

◆ AtRestGetPhysicalInteractionLength()

G4double G4eplusPolarizedAnnihilation::AtRestGetPhysicalInteractionLength ( const G4Track track,
G4ForceCondition condition 
)
inlinevirtual

Reimplemented from G4VDiscreteProcess.

Definition at line 135 of file G4eplusPolarizedAnnihilation.hh.

137{
139 return 0.0;
140}
G4double condition(const G4ErrorSymMatrix &m)
@ NotForced

◆ BuildAsymmetryTable()

void G4eplusPolarizedAnnihilation::BuildAsymmetryTable ( const G4ParticleDefinition part)

Definition at line 272 of file G4eplusPolarizedAnnihilation.cc.

273{
274 // Access to materials
275 const G4ProductionCutsTable* theCoupleTable=
277 size_t numOfCouples = theCoupleTable->GetTableSize();
278 G4cout<<" annih-numOfCouples="<<numOfCouples<<"\n";
279 for(size_t i=0; i<numOfCouples; ++i) {
280 G4cout<<"annih- "<<i<<"/"<<numOfCouples<<"\n";
281 if (!theAsymmetryTable) break;
282 G4cout<<"annih- "<<theAsymmetryTable->GetFlag(i)<<"\n";
283 if (theAsymmetryTable->GetFlag(i)) {
284 G4cout<<" building pol-annih ... \n";
285
286 // create physics vector and fill it
287 const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(i);
288
289 // use same parameters as for lambda
290 G4PhysicsVector* aVector = LambdaPhysicsVector(couple);
291 G4PhysicsVector* tVector = LambdaPhysicsVector(couple);
292
293 for (G4int j = 0 ; j < LambdaBinning() ; ++j ) {
294 G4double lowEdgeEnergy = aVector->GetLowEdgeEnergy(j);
295 G4double tasm=0.;
296 G4double asym = ComputeAsymmetry(lowEdgeEnergy, couple, part, 0., tasm);
297 aVector->PutValue(j,asym);
298 tVector->PutValue(j,tasm);
299 }
300
301 G4PhysicsTableHelper::SetPhysicsVector(theAsymmetryTable, i, aVector);
302 G4PhysicsTableHelper::SetPhysicsVector(theTransverseAsymmetryTable, i, tVector);
303 }
304 }
305
306}
int G4int
Definition: G4Types.hh:66
G4DLLIMPORT std::ostream G4cout
static void SetPhysicsVector(G4PhysicsTable *physTable, size_t idx, G4PhysicsVector *vec)
G4bool GetFlag(size_t i) const
virtual G4double GetLowEdgeEnergy(size_t binNumber) const
void PutValue(size_t index, G4double theValue)
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
static G4ProductionCutsTable * GetProductionCutsTable()
G4PhysicsVector * LambdaPhysicsVector(const G4MaterialCutsCouple *)
G4int LambdaBinning() const
G4double ComputeAsymmetry(G4double energy, const G4MaterialCutsCouple *couple, const G4ParticleDefinition &particle, G4double cut, G4double &tasm)

Referenced by BuildPhysicsTable().

◆ BuildPhysicsTable()

void G4eplusPolarizedAnnihilation::BuildPhysicsTable ( const G4ParticleDefinition pd)
virtual

Reimplemented from G4VProcess.

Definition at line 258 of file G4eplusPolarizedAnnihilation.cc.

259{
262}
void BuildPhysicsTable(const G4ParticleDefinition &)
void BuildAsymmetryTable(const G4ParticleDefinition &part)

◆ ComputeAsymmetry()

G4double G4eplusPolarizedAnnihilation::ComputeAsymmetry ( G4double  energy,
const G4MaterialCutsCouple couple,
const G4ParticleDefinition particle,
G4double  cut,
G4double tasm 
)

Definition at line 311 of file G4eplusPolarizedAnnihilation.cc.

316{
317 G4double lAsymmetry = 0.0;
318 tAsymmetry = 0.0;
319
320 // calculate polarized cross section
321 theTargetPolarization=G4ThreeVector(0.,0.,1.);
322 emModel->SetTargetPolarization(theTargetPolarization);
323 emModel->SetBeamPolarization(theTargetPolarization);
324 G4double sigma2=emModel->CrossSection(couple,&aParticle,energy,cut,energy);
325
326 // calculate transversely polarized cross section
327 theTargetPolarization=G4ThreeVector(1.,0.,0.);
328 emModel->SetTargetPolarization(theTargetPolarization);
329 emModel->SetBeamPolarization(theTargetPolarization);
330 G4double sigma3=emModel->CrossSection(couple,&aParticle,energy,cut,energy);
331
332 // calculate unpolarized cross section
333 theTargetPolarization=G4ThreeVector();
334 emModel->SetTargetPolarization(theTargetPolarization);
335 emModel->SetBeamPolarization(theTargetPolarization);
336 G4double sigma0=emModel->CrossSection(couple,&aParticle,energy,cut,energy);
337
338 // determine assymmetries
339 if (sigma0>0.) {
340 lAsymmetry=sigma2/sigma0-1.;
341 tAsymmetry=sigma3/sigma0-1.;
342 }
343 return lAsymmetry;
344
345}
CLHEP::Hep3Vector G4ThreeVector
void SetTargetPolarization(const G4ThreeVector &pTarget)
void SetBeamPolarization(const G4ThreeVector &pBeam)
G4double CrossSection(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:418

Referenced by BuildAsymmetryTable().

◆ GetMeanFreePath()

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

Implements G4VDiscreteProcess.

Definition at line 126 of file G4eplusPolarizedAnnihilation.cc.

129{
130 G4double mfp = G4VEmProcess::GetMeanFreePath(track, previousStepSize, condition);
131
132 if (theAsymmetryTable) {
133
134 G4Material* aMaterial = track.GetMaterial();
135 G4VPhysicalVolume* aPVolume = track.GetVolume();
136 G4LogicalVolume* aLVolume = aPVolume->GetLogicalVolume();
137
138 // G4Material* bMaterial = aLVolume->GetMaterial();
140
141 const G4bool volumeIsPolarized = polarizationManger->IsPolarized(aLVolume);
142 G4StokesVector electronPolarization = polarizationManger->GetVolumePolarization(aLVolume);
143
144 if (!volumeIsPolarized || mfp == DBL_MAX) return mfp;
145
146 // *** get asymmetry, if target is polarized ***
147 const G4DynamicParticle* aDynamicPositron = track.GetDynamicParticle();
148 const G4double positronEnergy = aDynamicPositron->GetKineticEnergy();
149 const G4StokesVector positronPolarization = track.GetPolarization();
150 const G4ParticleMomentum positronDirection0 = aDynamicPositron->GetMomentumDirection();
151
152 if (verboseLevel>=2) {
153
154 G4cout << " Mom " << positronDirection0 << G4endl;
155 G4cout << " Polarization " << positronPolarization << G4endl;
156 G4cout << " MaterialPol. " << electronPolarization << G4endl;
157 G4cout << " Phys. Volume " << aPVolume->GetName() << G4endl;
158 G4cout << " Log. Volume " << aLVolume->GetName() << G4endl;
159 G4cout << " Material " << aMaterial << G4endl;
160 }
161
162 G4bool isOutRange;
164 G4double lAsymmetry = (*theAsymmetryTable)(idx)->
165 GetValue(positronEnergy, isOutRange);
166 G4double tAsymmetry = (*theTransverseAsymmetryTable)(idx)->
167 GetValue(positronEnergy, isOutRange);
168
169 G4double polZZ = positronPolarization.z()*
170 electronPolarization*positronDirection0;
171 G4double polXX = positronPolarization.x()*
172 electronPolarization*G4PolarizationHelper::GetParticleFrameX(positronDirection0);
173 G4double polYY = positronPolarization.y()*
174 electronPolarization*G4PolarizationHelper::GetParticleFrameY(positronDirection0);
175
176 G4double impact = 1. + polZZ*lAsymmetry + (polXX + polYY)*tAsymmetry;
177
178 mfp *= 1. / impact;
179
180 if (verboseLevel>=2) {
181 G4cout << " MeanFreePath: " << mfp / mm << " mm " << G4endl;
182 G4cout << " Asymmetry: " << lAsymmetry << ", " << tAsymmetry << G4endl;
183 G4cout << " PolProduct: " << polXX << ", " << polYY << ", " << polZZ << G4endl;
184 }
185 }
186
187 return mfp;
188}
bool G4bool
Definition: G4Types.hh:67
#define G4endl
Definition: G4ios.hh:52
double z() const
double x() const
double y() const
const G4ThreeVector & GetMomentumDirection() const
G4double GetKineticEnergy() const
G4String GetName() const
static G4ThreeVector GetParticleFrameY(const G4ThreeVector &)
static G4ThreeVector GetParticleFrameX(const G4ThreeVector &)
bool IsPolarized(G4LogicalVolume *lVol) const
static G4PolarizationManager * GetInstance()
const G4ThreeVector & GetVolumePolarization(G4LogicalVolume *lVol) const
G4VPhysicalVolume * GetVolume() const
G4Material * GetMaterial() const
const G4DynamicParticle * GetDynamicParticle() const
const G4ThreeVector & GetPolarization() const
G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
size_t CurrentMaterialCutsCoupleIndex() const
G4LogicalVolume * GetLogicalVolume() const
const G4String & GetName() const
G4int verboseLevel
Definition: G4VProcess.hh:368
#define DBL_MAX
Definition: templates.hh:83

◆ InitialiseProcess()

void G4eplusPolarizedAnnihilation::InitialiseProcess ( const G4ParticleDefinition )
protectedvirtual

Implements G4VEmProcess.

Definition at line 101 of file G4eplusPolarizedAnnihilation.cc.

102{
103 if(!isInitialised) {
104 isInitialised = true;
105 // SetVerboseLevel(3);
106 SetBuildTableFlag(true);
109 G4double emin = 0.1*keV;
110 G4double emax = 100.*TeV;
111 SetLambdaBinning(120);
112 SetMinKinEnergy(emin);
113 SetMaxKinEnergy(emax);
114 emModel = new G4PolarizedAnnihilationModel();
115 emModel->SetLowEnergyLimit(emin);
116 emModel->SetHighEnergyLimit(emax);
117 AddEmModel(1, emModel);
118 }
119}
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:585
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:592
void AddEmModel(G4int, G4VEmModel *, const G4Region *region=0)
void SetBuildTableFlag(G4bool val)
void SetMinKinEnergy(G4double e)
void SetLambdaBinning(G4int nbins)
void SetSecondaryParticle(const G4ParticleDefinition *p)
void SetStartFromNullFlag(G4bool val)
void SetMaxKinEnergy(G4double e)

◆ IsApplicable()

G4bool G4eplusPolarizedAnnihilation::IsApplicable ( const G4ParticleDefinition p)
inlinevirtual

Implements G4VEmProcess.

Definition at line 128 of file G4eplusPolarizedAnnihilation.hh.

129{
130 return (&p == G4Positron::Positron());
131}
static G4Positron * Positron()
Definition: G4Positron.cc:94

◆ PostStepGetPhysicalInteractionLength()

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

Reimplemented from G4VDiscreteProcess.

Definition at line 192 of file G4eplusPolarizedAnnihilation.cc.

196{
198
199 if (theAsymmetryTable) {
200
201 G4Material* aMaterial = track.GetMaterial();
202 G4VPhysicalVolume* aPVolume = track.GetVolume();
203 G4LogicalVolume* aLVolume = aPVolume->GetLogicalVolume();
204
205 // G4Material* bMaterial = aLVolume->GetMaterial();
207
208 const G4bool volumeIsPolarized = polarizationManger->IsPolarized(aLVolume);
209 G4StokesVector electronPolarization = polarizationManger->GetVolumePolarization(aLVolume);
210
211 if (!volumeIsPolarized || mfp == DBL_MAX) return mfp;
212
213 // *** get asymmetry, if target is polarized ***
214 const G4DynamicParticle* aDynamicPositron = track.GetDynamicParticle();
215 const G4double positronEnergy = aDynamicPositron->GetKineticEnergy();
216 const G4StokesVector positronPolarization = track.GetPolarization();
217 const G4ParticleMomentum positronDirection0 = aDynamicPositron->GetMomentumDirection();
218
219 if (verboseLevel>=2) {
220
221 G4cout << " Mom " << positronDirection0 << G4endl;
222 G4cout << " Polarization " << positronPolarization << G4endl;
223 G4cout << " MaterialPol. " << electronPolarization << G4endl;
224 G4cout << " Phys. Volume " << aPVolume->GetName() << G4endl;
225 G4cout << " Log. Volume " << aLVolume->GetName() << G4endl;
226 G4cout << " Material " << aMaterial << G4endl;
227 }
228
229 G4bool isOutRange;
231 G4double lAsymmetry = (*theAsymmetryTable)(idx)->
232 GetValue(positronEnergy, isOutRange);
233 G4double tAsymmetry = (*theTransverseAsymmetryTable)(idx)->
234 GetValue(positronEnergy, isOutRange);
235
236 G4double polZZ = positronPolarization.z()*
237 electronPolarization*positronDirection0;
238 G4double polXX = positronPolarization.x()*
239 electronPolarization*G4PolarizationHelper::GetParticleFrameX(positronDirection0);
240 G4double polYY = positronPolarization.y()*
241 electronPolarization*G4PolarizationHelper::GetParticleFrameY(positronDirection0);
242
243 G4double impact = 1. + polZZ*lAsymmetry + (polXX + polYY)*tAsymmetry;
244
245 mfp *= 1. / impact;
246
247 if (verboseLevel>=2) {
248 G4cout << " MeanFreePath: " << mfp / mm << " mm " << G4endl;
249 G4cout << " Asymmetry: " << lAsymmetry << ", " << tAsymmetry << G4endl;
250 G4cout << " PolProduct: " << polXX << ", " << polYY << ", " << polZZ << G4endl;
251 }
252 }
253
254 return mfp;
255}
G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)

◆ PreparePhysicsTable()

void G4eplusPolarizedAnnihilation::PreparePhysicsTable ( const G4ParticleDefinition pd)
virtual

Reimplemented from G4VProcess.

Definition at line 265 of file G4eplusPolarizedAnnihilation.cc.

266{
268 theAsymmetryTable = G4PhysicsTableHelper::PreparePhysicsTable(theAsymmetryTable);
269 theTransverseAsymmetryTable = G4PhysicsTableHelper::PreparePhysicsTable(theTransverseAsymmetryTable);
270}
static G4PhysicsTable * PreparePhysicsTable(G4PhysicsTable *physTable)
void PreparePhysicsTable(const G4ParticleDefinition &)

◆ PrintInfo()

void G4eplusPolarizedAnnihilation::PrintInfo ( )
virtual

Implements G4VEmProcess.

Definition at line 351 of file G4eplusPolarizedAnnihilation.cc.

352{
353 G4cout << " Polarized model for annihilation into 2 photons"
354 << G4endl;
355}

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