G4eplusPolarizedAnnihilation Class Reference

#include <G4eplusPolarizedAnnihilation.hh>

Inheritance diagram for G4eplusPolarizedAnnihilation:

G4VEmProcess G4VDiscreteProcess G4VProcess

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)

Protected Member Functions

virtual void InitialiseProcess (const G4ParticleDefinition *)

Detailed Description

Definition at line 63 of file G4eplusPolarizedAnnihilation.hh.


Constructor & Destructor Documentation

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

Definition at line 75 of file G4eplusPolarizedAnnihilation.cc.

References G4VProcess::enableAtRestDoIt, fAnnihilation, and G4VProcess::SetProcessSubType().

00076   : G4VEmProcess(name), isInitialised(false),
00077     theAsymmetryTable(NULL),
00078     theTransverseAsymmetryTable(NULL)
00079 {
00080   enableAtRestDoIt = true;
00081   SetProcessSubType(fAnnihilation);
00082   emModel = 0; 
00083 }

G4eplusPolarizedAnnihilation::~G4eplusPolarizedAnnihilation (  )  [virtual]

Definition at line 87 of file G4eplusPolarizedAnnihilation.cc.

References G4PhysicsTable::clearAndDestroy().

00088 {
00089   if (theAsymmetryTable) {
00090     theAsymmetryTable->clearAndDestroy();
00091     delete theAsymmetryTable;
00092   }
00093   if (theTransverseAsymmetryTable) {
00094     theTransverseAsymmetryTable->clearAndDestroy();
00095     delete theTransverseAsymmetryTable;
00096   }
00097 }


Member Function Documentation

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

Reimplemented from G4VDiscreteProcess.

Definition at line 359 of file G4eplusPolarizedAnnihilation.cc.

References G4ParticleChangeForGamma::AddSecondary(), G4VEmProcess::fParticleChange, fStopAndKill, G4UniformRand, G4Gamma::Gamma(), G4ParticleChangeForGamma::InitializeForPostStep(), G4VParticleChange::ProposeTrackStatus(), and G4VParticleChange::SetNumberOfSecondaries().

00367         : Effects due to binding of atomic electrons are negliged.
00368 {
00369   fParticleChange.InitializeForPostStep(aTrack);
00370 
00371   fParticleChange.SetNumberOfSecondaries(2);
00372 
00373   G4double cosTeta = 2.*G4UniformRand()-1. , sinTeta = std::sqrt(1.-cosTeta*cosTeta);
00374   G4double phi     = twopi * G4UniformRand();
00375   G4ThreeVector direction (sinTeta*std::cos(phi), sinTeta*std::sin(phi), cosTeta);
00376   fParticleChange.AddSecondary( new G4DynamicParticle (G4Gamma::Gamma(),
00377                                             direction, electron_mass_c2) );
00378   fParticleChange.AddSecondary( new G4DynamicParticle (G4Gamma::Gamma(),
00379                                            -direction, electron_mass_c2) );
00380   // Kill the incident positron
00381   //
00382   fParticleChange.ProposeTrackStatus(fStopAndKill);
00383   return &fParticleChange;
00384 }

G4double G4eplusPolarizedAnnihilation::AtRestGetPhysicalInteractionLength ( const G4Track track,
G4ForceCondition condition 
) [inline, virtual]

Reimplemented from G4VDiscreteProcess.

Definition at line 135 of file G4eplusPolarizedAnnihilation.hh.

References NotForced.

00137 {
00138   *condition = NotForced;
00139   return 0.0;
00140 }

void G4eplusPolarizedAnnihilation::BuildAsymmetryTable ( const G4ParticleDefinition part  ) 

Definition at line 272 of file G4eplusPolarizedAnnihilation.cc.

References ComputeAsymmetry(), G4cout, G4PhysicsTable::GetFlag(), G4PhysicsVector::GetLowEdgeEnergy(), G4ProductionCutsTable::GetMaterialCutsCouple(), G4ProductionCutsTable::GetProductionCutsTable(), G4ProductionCutsTable::GetTableSize(), G4VEmProcess::LambdaBinning(), G4VEmProcess::LambdaPhysicsVector(), G4PhysicsVector::PutValue(), and G4PhysicsTableHelper::SetPhysicsVector().

Referenced by BuildPhysicsTable().

00273 {
00274   // Access to materials
00275   const G4ProductionCutsTable* theCoupleTable=
00276         G4ProductionCutsTable::GetProductionCutsTable();
00277   size_t numOfCouples = theCoupleTable->GetTableSize();
00278   G4cout<<" annih-numOfCouples="<<numOfCouples<<"\n";
00279   for(size_t i=0; i<numOfCouples; ++i) {
00280     G4cout<<"annih- "<<i<<"/"<<numOfCouples<<"\n";
00281     if (!theAsymmetryTable) break;
00282     G4cout<<"annih- "<<theAsymmetryTable->GetFlag(i)<<"\n";
00283     if (theAsymmetryTable->GetFlag(i)) {
00284      G4cout<<" building pol-annih ... \n";
00285 
00286       // create physics vector and fill it
00287       const G4MaterialCutsCouple* couple = theCoupleTable->GetMaterialCutsCouple(i);
00288 
00289       // use same parameters as for lambda
00290       G4PhysicsVector* aVector = LambdaPhysicsVector(couple);
00291       G4PhysicsVector* tVector = LambdaPhysicsVector(couple);
00292 
00293       for (G4int j = 0 ; j < LambdaBinning() ; ++j ) {
00294         G4double lowEdgeEnergy = aVector->GetLowEdgeEnergy(j);
00295         G4double tasm=0.;
00296         G4double asym = ComputeAsymmetry(lowEdgeEnergy, couple, part, 0., tasm);
00297         aVector->PutValue(j,asym);
00298         tVector->PutValue(j,tasm);
00299       }
00300 
00301       G4PhysicsTableHelper::SetPhysicsVector(theAsymmetryTable, i, aVector);
00302       G4PhysicsTableHelper::SetPhysicsVector(theTransverseAsymmetryTable, i, tVector);
00303     }
00304   }
00305 
00306 }

void G4eplusPolarizedAnnihilation::BuildPhysicsTable ( const G4ParticleDefinition  )  [virtual]

Reimplemented from G4VEmProcess.

Definition at line 258 of file G4eplusPolarizedAnnihilation.cc.

References BuildAsymmetryTable(), and G4VEmProcess::BuildPhysicsTable().

00259 {
00260   G4VEmProcess::BuildPhysicsTable(pd);
00261   BuildAsymmetryTable(pd);  
00262 }

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

Definition at line 311 of file G4eplusPolarizedAnnihilation.cc.

References G4VEmModel::CrossSection(), G4PolarizedAnnihilationModel::SetBeamPolarization(), and G4PolarizedAnnihilationModel::SetTargetPolarization().

Referenced by BuildAsymmetryTable().

00316 {
00317  G4double lAsymmetry = 0.0;
00318           tAsymmetry = 0.0;
00319 
00320  // calculate polarized cross section
00321  theTargetPolarization=G4ThreeVector(0.,0.,1.);
00322  emModel->SetTargetPolarization(theTargetPolarization);
00323  emModel->SetBeamPolarization(theTargetPolarization);
00324  G4double sigma2=emModel->CrossSection(couple,&aParticle,energy,cut,energy);
00325 
00326  // calculate transversely polarized cross section
00327  theTargetPolarization=G4ThreeVector(1.,0.,0.);
00328  emModel->SetTargetPolarization(theTargetPolarization);
00329  emModel->SetBeamPolarization(theTargetPolarization);
00330  G4double sigma3=emModel->CrossSection(couple,&aParticle,energy,cut,energy);
00331 
00332  // calculate unpolarized cross section
00333  theTargetPolarization=G4ThreeVector();
00334  emModel->SetTargetPolarization(theTargetPolarization);
00335  emModel->SetBeamPolarization(theTargetPolarization);
00336  G4double sigma0=emModel->CrossSection(couple,&aParticle,energy,cut,energy);
00337 
00338  // determine assymmetries
00339   if (sigma0>0.) {
00340     lAsymmetry=sigma2/sigma0-1.;
00341     tAsymmetry=sigma3/sigma0-1.;
00342                  }
00343  return lAsymmetry;
00344 
00345 }

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

Reimplemented from G4VEmProcess.

Definition at line 126 of file G4eplusPolarizedAnnihilation.cc.

References G4VEmProcess::CurrentMaterialCutsCoupleIndex(), DBL_MAX, G4cout, G4endl, G4Track::GetDynamicParticle(), G4PolarizationManager::GetInstance(), G4DynamicParticle::GetKineticEnergy(), G4VPhysicalVolume::GetLogicalVolume(), G4Track::GetMaterial(), G4VEmProcess::GetMeanFreePath(), G4DynamicParticle::GetMomentumDirection(), G4LogicalVolume::GetName(), G4VPhysicalVolume::GetName(), G4PolarizationHelper::GetParticleFrameX(), G4PolarizationHelper::GetParticleFrameY(), G4Track::GetPolarization(), G4Track::GetVolume(), G4PolarizationManager::GetVolumePolarization(), G4PolarizationManager::IsPolarized(), and G4VProcess::verboseLevel.

00129 {
00130   G4double mfp = G4VEmProcess::GetMeanFreePath(track, previousStepSize, condition);
00131 
00132   if (theAsymmetryTable) {
00133 
00134     G4Material*         aMaterial = track.GetMaterial();
00135     G4VPhysicalVolume*  aPVolume  = track.GetVolume();
00136     G4LogicalVolume*    aLVolume  = aPVolume->GetLogicalVolume();
00137     
00138     //   G4Material* bMaterial = aLVolume->GetMaterial();
00139     G4PolarizationManager * polarizationManger = G4PolarizationManager::GetInstance();
00140     
00141     const G4bool volumeIsPolarized = polarizationManger->IsPolarized(aLVolume);
00142     G4StokesVector electronPolarization = polarizationManger->GetVolumePolarization(aLVolume);
00143 
00144     if (!volumeIsPolarized || mfp == DBL_MAX) return mfp;
00145      
00146     // *** get asymmetry, if target is polarized ***
00147     const G4DynamicParticle* aDynamicPositron = track.GetDynamicParticle();
00148     const G4double positronEnergy = aDynamicPositron->GetKineticEnergy();
00149     const G4StokesVector positronPolarization = track.GetPolarization();
00150     const G4ParticleMomentum positronDirection0 = aDynamicPositron->GetMomentumDirection();
00151 
00152     if (verboseLevel>=2) {
00153       
00154       G4cout << " Mom " << positronDirection0  << G4endl;
00155       G4cout << " Polarization " << positronPolarization  << G4endl;
00156       G4cout << " MaterialPol. " << electronPolarization  << G4endl;
00157       G4cout << " Phys. Volume " << aPVolume->GetName() << G4endl;
00158       G4cout << " Log. Volume  " << aLVolume->GetName() << G4endl;
00159       G4cout << " Material     " << aMaterial          << G4endl;
00160     }
00161     
00162     G4bool isOutRange;
00163     G4int idx= CurrentMaterialCutsCoupleIndex();
00164     G4double lAsymmetry = (*theAsymmetryTable)(idx)->
00165                                   GetValue(positronEnergy, isOutRange);
00166     G4double tAsymmetry = (*theTransverseAsymmetryTable)(idx)->
00167                                   GetValue(positronEnergy, isOutRange);
00168 
00169     G4double polZZ = positronPolarization.z()*
00170       electronPolarization*positronDirection0;
00171     G4double polXX = positronPolarization.x()*
00172       electronPolarization*G4PolarizationHelper::GetParticleFrameX(positronDirection0);
00173     G4double polYY = positronPolarization.y()*
00174       electronPolarization*G4PolarizationHelper::GetParticleFrameY(positronDirection0);
00175 
00176     G4double impact = 1. + polZZ*lAsymmetry + (polXX + polYY)*tAsymmetry;
00177 
00178     mfp *= 1. / impact;
00179 
00180     if (verboseLevel>=2) {
00181       G4cout << " MeanFreePath:  " << mfp / mm << " mm " << G4endl;
00182       G4cout << " Asymmetry:     " << lAsymmetry << ", " << tAsymmetry  << G4endl;
00183       G4cout << " PolProduct:    " << polXX << ", " << polYY << ", " << polZZ << G4endl;
00184     }
00185   }
00186   
00187   return mfp;
00188 }

void G4eplusPolarizedAnnihilation::InitialiseProcess ( const G4ParticleDefinition  )  [protected, virtual]

Implements G4VEmProcess.

Definition at line 101 of file G4eplusPolarizedAnnihilation.cc.

References G4VEmProcess::AddEmModel(), G4Gamma::Gamma(), G4VEmProcess::SetBuildTableFlag(), G4VEmModel::SetHighEnergyLimit(), G4VEmProcess::SetLambdaBinning(), G4VEmModel::SetLowEnergyLimit(), G4VEmProcess::SetMaxKinEnergy(), G4VEmProcess::SetMinKinEnergy(), G4VEmProcess::SetSecondaryParticle(), and G4VEmProcess::SetStartFromNullFlag().

00102 {
00103   if(!isInitialised) {
00104     isInitialised = true;
00105     //    SetVerboseLevel(3);
00106     SetBuildTableFlag(true);
00107     SetStartFromNullFlag(false);
00108     SetSecondaryParticle(G4Gamma::Gamma());
00109     G4double emin = 0.1*keV;
00110     G4double emax = 100.*TeV;
00111     SetLambdaBinning(120);
00112     SetMinKinEnergy(emin);
00113     SetMaxKinEnergy(emax);
00114     emModel = new G4PolarizedAnnihilationModel();
00115     emModel->SetLowEnergyLimit(emin);
00116     emModel->SetHighEnergyLimit(emax);
00117     AddEmModel(1, emModel);
00118   }
00119 }

G4bool G4eplusPolarizedAnnihilation::IsApplicable ( const G4ParticleDefinition p  )  [inline, virtual]

Implements G4VEmProcess.

Definition at line 128 of file G4eplusPolarizedAnnihilation.hh.

References G4Positron::Positron().

00129 {
00130   return (&p == G4Positron::Positron());
00131 }

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

Reimplemented from G4VEmProcess.

Definition at line 192 of file G4eplusPolarizedAnnihilation.cc.

References G4VEmProcess::CurrentMaterialCutsCoupleIndex(), DBL_MAX, G4cout, G4endl, G4Track::GetDynamicParticle(), G4PolarizationManager::GetInstance(), G4DynamicParticle::GetKineticEnergy(), G4VPhysicalVolume::GetLogicalVolume(), G4Track::GetMaterial(), G4DynamicParticle::GetMomentumDirection(), G4LogicalVolume::GetName(), G4VPhysicalVolume::GetName(), G4PolarizationHelper::GetParticleFrameX(), G4PolarizationHelper::GetParticleFrameY(), G4Track::GetPolarization(), G4Track::GetVolume(), G4PolarizationManager::GetVolumePolarization(), G4PolarizationManager::IsPolarized(), G4VEmProcess::PostStepGetPhysicalInteractionLength(), and G4VProcess::verboseLevel.

00196 {
00197   G4double mfp = G4VEmProcess::PostStepGetPhysicalInteractionLength(track, previousStepSize, condition);
00198 
00199   if (theAsymmetryTable) {
00200 
00201     G4Material*         aMaterial = track.GetMaterial();
00202     G4VPhysicalVolume*  aPVolume  = track.GetVolume();
00203     G4LogicalVolume*    aLVolume  = aPVolume->GetLogicalVolume();
00204     
00205     //   G4Material* bMaterial = aLVolume->GetMaterial();
00206     G4PolarizationManager * polarizationManger = G4PolarizationManager::GetInstance();
00207     
00208     const G4bool volumeIsPolarized = polarizationManger->IsPolarized(aLVolume);
00209     G4StokesVector electronPolarization = polarizationManger->GetVolumePolarization(aLVolume);
00210 
00211     if (!volumeIsPolarized || mfp == DBL_MAX) return mfp;
00212      
00213     // *** get asymmetry, if target is polarized ***
00214     const G4DynamicParticle* aDynamicPositron = track.GetDynamicParticle();
00215     const G4double positronEnergy = aDynamicPositron->GetKineticEnergy();
00216     const G4StokesVector positronPolarization = track.GetPolarization();
00217     const G4ParticleMomentum positronDirection0 = aDynamicPositron->GetMomentumDirection();
00218 
00219     if (verboseLevel>=2) {
00220       
00221       G4cout << " Mom " << positronDirection0  << G4endl;
00222       G4cout << " Polarization " << positronPolarization  << G4endl;
00223       G4cout << " MaterialPol. " << electronPolarization  << G4endl;
00224       G4cout << " Phys. Volume " << aPVolume->GetName() << G4endl;
00225       G4cout << " Log. Volume  " << aLVolume->GetName() << G4endl;
00226       G4cout << " Material     " << aMaterial          << G4endl;
00227     }
00228     
00229     G4bool isOutRange;
00230     G4int idx= CurrentMaterialCutsCoupleIndex();
00231     G4double lAsymmetry = (*theAsymmetryTable)(idx)->
00232                                   GetValue(positronEnergy, isOutRange);
00233     G4double tAsymmetry = (*theTransverseAsymmetryTable)(idx)->
00234                                   GetValue(positronEnergy, isOutRange);
00235 
00236     G4double polZZ = positronPolarization.z()*
00237       electronPolarization*positronDirection0;
00238     G4double polXX = positronPolarization.x()*
00239       electronPolarization*G4PolarizationHelper::GetParticleFrameX(positronDirection0);
00240     G4double polYY = positronPolarization.y()*
00241       electronPolarization*G4PolarizationHelper::GetParticleFrameY(positronDirection0);
00242 
00243     G4double impact = 1. + polZZ*lAsymmetry + (polXX + polYY)*tAsymmetry;
00244 
00245     mfp *= 1. / impact;
00246 
00247     if (verboseLevel>=2) {
00248       G4cout << " MeanFreePath:  " << mfp / mm << " mm " << G4endl;
00249       G4cout << " Asymmetry:     " << lAsymmetry << ", " << tAsymmetry  << G4endl;
00250       G4cout << " PolProduct:    " << polXX << ", " << polYY << ", " << polZZ << G4endl;
00251     }
00252   }
00253   
00254   return mfp;
00255 }

void G4eplusPolarizedAnnihilation::PreparePhysicsTable ( const G4ParticleDefinition  )  [virtual]

Reimplemented from G4VEmProcess.

Definition at line 265 of file G4eplusPolarizedAnnihilation.cc.

References G4PhysicsTableHelper::PreparePhysicsTable(), and G4VEmProcess::PreparePhysicsTable().

00266 {
00267   G4VEmProcess::PreparePhysicsTable(pd);
00268   theAsymmetryTable = G4PhysicsTableHelper::PreparePhysicsTable(theAsymmetryTable);
00269   theTransverseAsymmetryTable = G4PhysicsTableHelper::PreparePhysicsTable(theTransverseAsymmetryTable);
00270 }

void G4eplusPolarizedAnnihilation::PrintInfo (  )  [virtual]

Implements G4VEmProcess.

Definition at line 351 of file G4eplusPolarizedAnnihilation.cc.

References G4cout, and G4endl.

00352 {
00353   G4cout << "      Polarized model for annihilation into 2 photons"
00354          << G4endl;
00355 }


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:51:54 2013 for Geant4 by  doxygen 1.4.7