Geant4  10.02.p03
G4eplusPolarizedAnnihilation Class Reference

#include <G4eplusPolarizedAnnihilation.hh>

Inheritance diagram for G4eplusPolarizedAnnihilation:
Collaboration diagram for G4eplusPolarizedAnnihilation:

Public Member Functions

 G4eplusPolarizedAnnihilation (const G4String &name="pol-annihil")
 
virtual ~G4eplusPolarizedAnnihilation ()
 
virtual void PrintInfo ()
 
virtual G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual G4double PostStepGetPhysicalInteractionLength (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
 
virtual void BuildPhysicsTable (const G4ParticleDefinition &)
 
- Public Member Functions inherited from G4eplusAnnihilation
 G4eplusAnnihilation (const G4String &name="annihil")
 
virtual ~G4eplusAnnihilation ()
 
virtual G4bool IsApplicable (const G4ParticleDefinition &p)
 
virtual G4VParticleChange * AtRestDoIt (const G4Track &track, const G4Step &stepData)
 
virtual G4double AtRestGetPhysicalInteractionLength (const G4Track &track, G4ForceCondition *condition)
 
- Public Member Functions inherited from G4VEmProcess
 G4VEmProcess (const G4String &name, G4ProcessType type=fElectromagnetic)
 
virtual ~G4VEmProcess ()
 
virtual void ProcessDescription (std::ostream &outFile) const
 
virtual void PreparePhysicsTable (const G4ParticleDefinition &)
 
virtual void StartTracking (G4Track *)
 
virtual G4VParticleChange * PostStepDoIt (const G4Track &, const G4Step &)
 
virtual G4bool StorePhysicsTable (const G4ParticleDefinition *, const G4String &directory, G4bool ascii=false)
 
virtual 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)
 
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=0)
 
G4VEmModelEmModel (G4int index=1) const
 
void SetEmModel (G4VEmModel *, G4int index=1)
 
void UpdateEmModel (const G4String &, G4double, G4double)
 
G4VEmModelGetModelByIndex (G4int idx=0, G4bool ver=false) const
 
const G4ElementGetCurrentElement () const
 
const G4VEmModelGetCurrentModel () 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 SetIntegral (G4bool val)
 
void SetBuildTableFlag (G4bool val)
 
- Public Member Functions inherited from G4VDiscreteProcess
 G4VDiscreteProcess (const G4String &, G4ProcessType aType=fNotDefined)
 
 G4VDiscreteProcess (G4VDiscreteProcess &)
 
virtual ~G4VDiscreteProcess ()
 
virtual G4double AlongStepGetPhysicalInteractionLength (const G4Track &, G4double, G4double, G4double &, G4GPILSelection *)
 
virtual G4VParticleChange * AlongStepDoIt (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
 
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 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 &)
 

Private Member Functions

void CleanTables ()
 
void BuildAsymmetryTables (const G4ParticleDefinition &part)
 
G4double ComputeAsymmetry (G4double energy, const G4MaterialCutsCouple *couple, const G4ParticleDefinition &particle, G4double cut, G4double &tasm)
 
G4double ComputeSaturationFactor (const G4Track &aTrack)
 
G4eplusPolarizedAnnihilationoperator= (const G4eplusPolarizedAnnihilation &right)
 
 G4eplusPolarizedAnnihilation (const G4eplusPolarizedAnnihilation &)
 

Private Attributes

G4bool isInitialised
 
G4PolarizedAnnihilationModelemModel
 
G4ThreeVector theTargetPolarization
 
G4PhysicsTabletheAsymmetryTable
 
G4PhysicsTabletheTransverseAsymmetryTable
 

Additional Inherited Members

- Static Public Member Functions inherited from G4VProcess
static const G4StringGetProcessTypeName (G4ProcessType)
 
- Protected Member Functions inherited from G4eplusAnnihilation
virtual void InitialiseProcess (const G4ParticleDefinition *)
 
- Protected Member Functions inherited from G4VEmProcess
virtual G4double MinPrimaryEnergy (const G4ParticleDefinition *, const G4Material *)
 
G4VEmModelSelectModel (G4double &kinEnergy, size_t index)
 
G4PhysicsVectorLambdaPhysicsVector (const G4MaterialCutsCouple *)
 
G4int LambdaBinning () const
 
G4double MinKinEnergy () const
 
G4double MaxKinEnergy () const
 
G4double PolarAngleLimit () const
 
G4bool IsIntegral () const
 
G4double RecalculateLambda (G4double kinEnergy, const G4MaterialCutsCouple *couple)
 
G4ParticleChangeForGamma * GetParticleChange ()
 
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)
 
const G4ElementGetTargetElement () const
 
const G4IsotopeGetTargetIsotope () const
 
- Protected Member Functions inherited from G4VProcess
void SubtractNumberOfInteractionLengthLeft (G4double previousStepSize)
 
void ClearNumberOfInteractionLengthLeft ()
 
- Protected Attributes inherited from G4VEmProcess
G4ParticleChangeForGamma fParticleChange
 
- Protected Attributes inherited from G4VProcess
const G4ProcessManageraProcessManager
 
G4VParticleChange * pParticleChange
 
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() [1/2]

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

Definition at line 74 of file G4eplusPolarizedAnnihilation.cc.

75  : G4eplusAnnihilation(name), isInitialised(false),
76  theAsymmetryTable(nullptr),
78 {
80  SetEmModel(emModel, 1);
81 }
void SetEmModel(G4VEmModel *, G4int index=1)
G4eplusAnnihilation(const G4String &name="annihil")
G4PolarizedAnnihilationModel * emModel
Here is the call graph for this function:

◆ ~G4eplusPolarizedAnnihilation()

G4eplusPolarizedAnnihilation::~G4eplusPolarizedAnnihilation ( )
virtual

Definition at line 85 of file G4eplusPolarizedAnnihilation.cc.

Here is the call graph for this function:

◆ G4eplusPolarizedAnnihilation() [2/2]

G4eplusPolarizedAnnihilation::G4eplusPolarizedAnnihilation ( const G4eplusPolarizedAnnihilation )
private

Member Function Documentation

◆ BuildAsymmetryTables()

void G4eplusPolarizedAnnihilation::BuildAsymmetryTables ( const G4ParticleDefinition part)
private

Definition at line 240 of file G4eplusPolarizedAnnihilation.cc.

242 {
243  // cleanup old, initialise new table
244  CleanTables();
249 
250  // Access to materials
251  const G4ProductionCutsTable* theCoupleTable=
253  size_t numOfCouples = theCoupleTable->GetTableSize();
254  //G4cout<<" annih-numOfCouples="<<numOfCouples<<"\n";
255  for(size_t i=0; i<numOfCouples; ++i) {
256  //G4cout<<"annih- "<<i<<"/"<<numOfCouples<<"\n";
257  if (!theAsymmetryTable) break;
258  //G4cout<<"annih- "<<theAsymmetryTable->GetFlag(i)<<"\n";
259  if (theAsymmetryTable->GetFlag(i)) {
260  //G4cout<<" building pol-annih ... \n";
261 
262  // create physics vector and fill it
263  const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(i);
264 
265  // use same parameters as for lambda
266  G4PhysicsVector* aVector = LambdaPhysicsVector(couple);
267  G4PhysicsVector* tVector = LambdaPhysicsVector(couple);
268 
269  for (G4int j = 0 ; j < LambdaBinning() ; ++j ) {
270  G4double lowEdgeEnergy = aVector->GetLowEdgeEnergy(j);
271  G4double tasm=0.;
272  G4double asym = ComputeAsymmetry(lowEdgeEnergy, couple, part, 0., tasm);
273  aVector->PutValue(j,asym);
274  tVector->PutValue(j,tasm);
275  }
278  }
279  }
280 }
static void SetPhysicsVector(G4PhysicsTable *physTable, size_t idx, G4PhysicsVector *vec)
G4int LambdaBinning() const
int G4int
Definition: G4Types.hh:78
static G4PhysicsTable * PreparePhysicsTable(G4PhysicsTable *physTable)
G4double ComputeAsymmetry(G4double energy, const G4MaterialCutsCouple *couple, const G4ParticleDefinition &particle, G4double cut, G4double &tasm)
G4double GetLowEdgeEnergy(size_t binNumber) const
void PutValue(size_t index, G4double theValue)
static G4ProductionCutsTable * GetProductionCutsTable()
G4PhysicsVector * LambdaPhysicsVector(const G4MaterialCutsCouple *)
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
G4bool GetFlag(size_t i) const
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ BuildPhysicsTable()

void G4eplusPolarizedAnnihilation::BuildPhysicsTable ( const G4ParticleDefinition part)
virtual

Reimplemented from G4VEmProcess.

Definition at line 227 of file G4eplusPolarizedAnnihilation.cc.

229 {
231  G4bool isMaster = true;
232  const G4eplusPolarizedAnnihilation* masterProcess =
233  static_cast<const G4eplusPolarizedAnnihilation*>(GetMasterProcess());
234  if(masterProcess && masterProcess != this) { isMaster = false; }
235  if(isMaster) { BuildAsymmetryTables(part); }
236 }
const G4VProcess * GetMasterProcess() const
Definition: G4VProcess.hh:538
virtual void BuildPhysicsTable(const G4ParticleDefinition &)
bool G4bool
Definition: G4Types.hh:79
void BuildAsymmetryTables(const G4ParticleDefinition &part)
Here is the call graph for this function:

◆ CleanTables()

void G4eplusPolarizedAnnihilation::CleanTables ( )
private

Definition at line 92 of file G4eplusPolarizedAnnihilation.cc.

Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeAsymmetry()

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

Definition at line 284 of file G4eplusPolarizedAnnihilation.cc.

289 {
290  G4double lAsymmetry = 0.0;
291  tAsymmetry = 0.0;
292 
293  // calculate polarized cross section
297  G4double sigma2=emModel->CrossSection(couple,&aParticle,energy,cut,energy);
298 
299  // calculate transversely polarized cross section
303  G4double sigma3=emModel->CrossSection(couple,&aParticle,energy,cut,energy);
304 
305  // calculate unpolarized cross section
309  G4double sigma0=emModel->CrossSection(couple,&aParticle,energy,cut,energy);
310 
311  // determine assymmetries
312  if (sigma0>0.) {
313  lAsymmetry=sigma2/sigma0-1.;
314  tAsymmetry=sigma3/sigma0-1.;
315  }
316  return lAsymmetry;
317 }
void SetBeamPolarization(const G4ThreeVector &pBeam)
CLHEP::Hep3Vector G4ThreeVector
double energy
Definition: plottest35.C:25
void SetTargetPolarization(const G4ThreeVector &pTarget)
G4double CrossSection(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:501
G4PolarizedAnnihilationModel * emModel
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeSaturationFactor()

G4double G4eplusPolarizedAnnihilation::ComputeSaturationFactor ( const G4Track &  aTrack)
private

Definition at line 157 of file G4eplusPolarizedAnnihilation.cc.

158 {
159  G4Material* aMaterial = track.GetMaterial();
160  G4VPhysicalVolume* aPVolume = track.GetVolume();
161  G4LogicalVolume* aLVolume = aPVolume->GetLogicalVolume();
162 
164 
165  const G4bool volumeIsPolarized = polarizationManger->IsPolarized(aLVolume);
166  G4StokesVector electronPolarization = polarizationManger->GetVolumePolarization(aLVolume);
167 
168  G4double factor = 1.0;
169 
170  if (volumeIsPolarized) {
171 
172  // *** get asymmetry, if target is polarized ***
173  const G4DynamicParticle* aDynamicPositron = track.GetDynamicParticle();
174  const G4double positronEnergy = aDynamicPositron->GetKineticEnergy();
175  const G4StokesVector positronPolarization = track.GetPolarization();
176  const G4ParticleMomentum positronDirection0 = aDynamicPositron->GetMomentumDirection();
177 
178  if (verboseLevel>=2) {
179  G4cout << "G4eplusPolarizedAnnihilation::ComputeSaturationFactor: " << G4endl;
180  G4cout << " Mom " << positronDirection0 << G4endl;
181  G4cout << " Polarization " << positronPolarization << G4endl;
182  G4cout << " MaterialPol. " << electronPolarization << G4endl;
183  G4cout << " Phys. Volume " << aPVolume->GetName() << G4endl;
184  G4cout << " Log. Volume " << aLVolume->GetName() << G4endl;
185  G4cout << " Material " << aMaterial << G4endl;
186  }
187 
188  size_t midx = CurrentMaterialCutsCoupleIndex();
189  const G4PhysicsVector* aVector = nullptr;
190  const G4PhysicsVector* bVector = nullptr;
191  if(midx < theAsymmetryTable->size()) {
192  aVector = (*theAsymmetryTable)(midx);
193  }
194  if(midx < theTransverseAsymmetryTable->size()) {
195  bVector = (*theTransverseAsymmetryTable)(midx);
196  }
197  if(aVector && bVector) {
198  G4double lAsymmetry = aVector->Value(positronEnergy);
199  G4double tAsymmetry = bVector->Value(positronEnergy);
200  G4double polZZ = positronPolarization.z()*
201  (electronPolarization*positronDirection0);
202  G4double polXX = positronPolarization.x()*
203  (electronPolarization*G4PolarizationHelper::GetParticleFrameX(positronDirection0));
204  G4double polYY = positronPolarization.y()*
205  (electronPolarization*G4PolarizationHelper::GetParticleFrameY(positronDirection0));
206 
207  factor /= (1. + polZZ*lAsymmetry + (polXX + polYY)*tAsymmetry);
208 
209  if (verboseLevel>=2) {
210  G4cout << " Asymmetry: " << lAsymmetry << ", " << tAsymmetry << G4endl;
211  G4cout << " PolProduct: " << polXX << ", " << polYY << ", " << polZZ << G4endl;
212  G4cout << " Factor: " << factor << G4endl;
213  }
214  } else {
216  ed << "Problem with asymmetry tables: material index " << midx
217  << " is out of range or tables are not filled";
218  G4Exception("G4eplusPolarizedAnnihilation::ComputeSaturationFactor","em0048",
219  JustWarning, ed, "");
220  }
221  }
222  return factor;
223 }
G4int verboseLevel
Definition: G4VProcess.hh:368
std::ostringstream G4ExceptionDescription
Definition: globals.hh:76
const G4ThreeVector & GetVolumePolarization(G4LogicalVolume *lVol) const
static G4PolarizationManager * GetInstance()
size_t CurrentMaterialCutsCoupleIndex() const
G4double GetKineticEnergy() const
bool IsPolarized(G4LogicalVolume *lVol) const
G4GLOB_DLL std::ostream G4cout
const G4String & GetName() const
bool G4bool
Definition: G4Types.hh:79
G4double Value(G4double theEnergy, size_t &lastidx) const
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
double x() const
const G4String & GetName() const
static const G4double factor
double y() const
const G4ThreeVector & GetMomentumDirection() const
double z() const
static G4ThreeVector GetParticleFrameY(const G4ThreeVector &)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4LogicalVolume * GetLogicalVolume() const
static G4ThreeVector GetParticleFrameX(const G4ThreeVector &)
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetMeanFreePath()

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

Reimplemented from G4VEmProcess.

Definition at line 108 of file G4eplusPolarizedAnnihilation.cc.

111 {
112  G4double mfp = G4VEmProcess::GetMeanFreePath(track, previousStepSize, condition);
113 
116  }
117  if (verboseLevel>=2) {
118  G4cout << "G4eplusPolarizedAnnihilation::MeanFreePath: "
119  << mfp / mm << " mm " << G4endl;
120  }
121  return mfp;
122 }
G4double condition(const G4ErrorSymMatrix &m)
G4int verboseLevel
Definition: G4VProcess.hh:368
G4double ComputeSaturationFactor(const G4Track &aTrack)
G4GLOB_DLL std::ostream G4cout
virtual G4double GetMeanFreePath(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
static const double mm
Definition: G4SIunits.hh:114
Here is the call graph for this function:

◆ operator=()

G4eplusPolarizedAnnihilation& G4eplusPolarizedAnnihilation::operator= ( const G4eplusPolarizedAnnihilation right)
private

◆ PostStepGetPhysicalInteractionLength()

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

Reimplemented from G4VEmProcess.

Definition at line 126 of file G4eplusPolarizedAnnihilation.cc.

130 {
131  // save previous value
133 
134  // *** compute uppolarized step limit ***
136  previousStepSize,
137  condition);
138 
141  if(nLength > 0.0) {
143  std::max(nLength - previousStepSize/curLength, 0.0);
144  }
145  x = theNumberOfInteractionLengthLeft * curLength;
146  }
147  if (verboseLevel>=2) {
148  G4cout << "G4eplusPolarizedAnnihilation::PostStepGetPhysicalInteractionLength: "
149  << x/mm << " mm " << G4endl;
150  }
151  return x;
152 }
G4double condition(const G4ErrorSymMatrix &m)
G4int verboseLevel
Definition: G4VProcess.hh:368
G4double ComputeSaturationFactor(const G4Track &aTrack)
G4double theNumberOfInteractionLengthLeft
Definition: G4VProcess.hh:293
virtual G4double PostStepGetPhysicalInteractionLength(const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)
G4GLOB_DLL std::ostream G4cout
G4double currentInteractionLength
Definition: G4VProcess.hh:297
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
static const double mm
Definition: G4SIunits.hh:114
Here is the call graph for this function:

◆ PrintInfo()

void G4eplusPolarizedAnnihilation::PrintInfo ( )
virtual

Reimplemented from G4eplusAnnihilation.

Definition at line 321 of file G4eplusPolarizedAnnihilation.cc.

322 {
323  G4cout << " Polarized model for annihilation into 2 photons"
324  << G4endl;
325 }
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61

Member Data Documentation

◆ emModel

G4PolarizedAnnihilationModel* G4eplusPolarizedAnnihilation::emModel
private

Definition at line 106 of file G4eplusPolarizedAnnihilation.hh.

◆ isInitialised

G4bool G4eplusPolarizedAnnihilation::isInitialised
private

Definition at line 103 of file G4eplusPolarizedAnnihilation.hh.

◆ theAsymmetryTable

G4PhysicsTable* G4eplusPolarizedAnnihilation::theAsymmetryTable
private

Definition at line 109 of file G4eplusPolarizedAnnihilation.hh.

◆ theTargetPolarization

G4ThreeVector G4eplusPolarizedAnnihilation::theTargetPolarization
private

Definition at line 107 of file G4eplusPolarizedAnnihilation.hh.

◆ theTransverseAsymmetryTable

G4PhysicsTable* G4eplusPolarizedAnnihilation::theTransverseAsymmetryTable
private

Definition at line 110 of file G4eplusPolarizedAnnihilation.hh.


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