G4eplusPolarizedAnnihilation.cc

Go to the documentation of this file.
00001 //
00002 // ********************************************************************
00003 // * License and Disclaimer                                           *
00004 // *                                                                  *
00005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
00006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
00007 // * conditions of the Geant4 Software License,  included in the file *
00008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
00009 // * include a list of copyright holders.                             *
00010 // *                                                                  *
00011 // * Neither the authors of this software system, nor their employing *
00012 // * institutes,nor the agencies providing financial support for this *
00013 // * work  make  any representation or  warranty, express or implied, *
00014 // * regarding  this  software system or assume any liability for its *
00015 // * use.  Please see the license in the file  LICENSE  and URL above *
00016 // * for the full disclaimer and the limitation of liability.         *
00017 // *                                                                  *
00018 // * This  code  implementation is the result of  the  scientific and *
00019 // * technical work of the GEANT4 collaboration.                      *
00020 // * By using,  copying,  modifying or  distributing the software (or *
00021 // * any work based  on the software)  you  agree  to acknowledge its *
00022 // * use  in  resulting  scientific  publications,  and indicate your *
00023 // * acceptance of all terms of the Geant4 Software license.          *
00024 // ********************************************************************
00025 //
00026 // $Id: G4eplusPolarizedAnnihilation.cc 69847 2013-05-16 09:36:18Z gcosmo $
00027 //
00028 // -------------------------------------------------------------------
00029 //
00030 // GEANT4 Class file
00031 //
00032 //
00033 // File name:     G4eplusPolarizedAnnihilation
00034 //
00035 // Author:        A. Schaelicke on base of Vladimir Ivanchenko / Michel Maire code
00036 //
00037 // Creation date: 02.07.2006
00038 //
00039 // Modifications:
00040 // 26-07-06 modified cross section  (P. Starovoitov)
00041 // 21-08-06 interface updated   (A. Schaelicke)
00042 // 11-06-07, add PostStepGetPhysicalInteractionLength (A.Schalicke)
00043 // 02-10-07, enable AtRest (V.Ivanchenko)
00044 //
00045 //
00046 // Class Description:
00047 //
00048 // Polarized process of e+ annihilation into 2 gammas
00049 //
00050 
00051 //
00052 // -------------------------------------------------------------------
00053 //
00054 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00055 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00056 
00057 #include "G4eplusPolarizedAnnihilation.hh"
00058 #include "G4PhysicalConstants.hh"
00059 #include "G4SystemOfUnits.hh"
00060 #include "G4MaterialCutsCouple.hh"
00061 #include "G4Gamma.hh"
00062 #include "G4PhysicsVector.hh"
00063 #include "G4PhysicsLogVector.hh"
00064 
00065 
00066 #include "G4PolarizedAnnihilationModel.hh"
00067 #include "G4PhysicsTableHelper.hh"
00068 #include "G4ProductionCutsTable.hh"
00069 #include "G4PolarizationManager.hh"
00070 #include "G4PolarizationHelper.hh"
00071 #include "G4StokesVector.hh"
00072 
00073 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00074 
00075 G4eplusPolarizedAnnihilation::G4eplusPolarizedAnnihilation(const G4String& name)
00076   : G4VEmProcess(name), isInitialised(false),
00077     theAsymmetryTable(NULL),
00078     theTransverseAsymmetryTable(NULL)
00079 {
00080   enableAtRestDoIt = true;
00081   SetProcessSubType(fAnnihilation);
00082   emModel = 0; 
00083 }
00084 
00085 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00086 
00087 G4eplusPolarizedAnnihilation::~G4eplusPolarizedAnnihilation()
00088 {
00089   if (theAsymmetryTable) {
00090     theAsymmetryTable->clearAndDestroy();
00091     delete theAsymmetryTable;
00092   }
00093   if (theTransverseAsymmetryTable) {
00094     theTransverseAsymmetryTable->clearAndDestroy();
00095     delete theTransverseAsymmetryTable;
00096   }
00097 }
00098 
00099 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00100 
00101 void G4eplusPolarizedAnnihilation::InitialiseProcess(const G4ParticleDefinition*)
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 }
00120 
00121 
00122 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00123 
00124   // for polarization
00125 
00126 G4double G4eplusPolarizedAnnihilation::GetMeanFreePath(const G4Track& track,
00127                               G4double previousStepSize,
00128                               G4ForceCondition* condition)
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 }
00189 
00190 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00191 
00192 G4double G4eplusPolarizedAnnihilation::PostStepGetPhysicalInteractionLength(
00193                               const G4Track& track,
00194                               G4double previousStepSize,
00195                               G4ForceCondition* condition)
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 }
00256 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00257 
00258 void G4eplusPolarizedAnnihilation::BuildPhysicsTable(const G4ParticleDefinition& pd) 
00259 {
00260   G4VEmProcess::BuildPhysicsTable(pd);
00261   BuildAsymmetryTable(pd);  
00262 }
00263 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00264 
00265 void G4eplusPolarizedAnnihilation::PreparePhysicsTable(const G4ParticleDefinition& pd)
00266 {
00267   G4VEmProcess::PreparePhysicsTable(pd);
00268   theAsymmetryTable = G4PhysicsTableHelper::PreparePhysicsTable(theAsymmetryTable);
00269   theTransverseAsymmetryTable = G4PhysicsTableHelper::PreparePhysicsTable(theTransverseAsymmetryTable);
00270 }
00271 
00272 void G4eplusPolarizedAnnihilation::BuildAsymmetryTable(const G4ParticleDefinition& part)
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 }
00307 
00308 
00309 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00310 
00311 G4double G4eplusPolarizedAnnihilation::ComputeAsymmetry(G4double energy,
00312                             const G4MaterialCutsCouple* couple,
00313                             const G4ParticleDefinition& aParticle,
00314                             G4double cut,
00315                             G4double &tAsymmetry)
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 }
00346 
00347 
00348 
00349 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00350 
00351 void G4eplusPolarizedAnnihilation::PrintInfo()
00352 {
00353   G4cout << "      Polarized model for annihilation into 2 photons"
00354          << G4endl;
00355 }
00356 
00357 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00358 
00359 G4VParticleChange* G4eplusPolarizedAnnihilation::AtRestDoIt(const G4Track& aTrack,
00360                                                      const G4Step& )
00361 //
00362 // Performs the e+ e- annihilation when both particles are assumed at rest.
00363 // It generates two back to back photons with energy = electron_mass.
00364 // The angular distribution is isotropic.
00365 // GEANT4 internal units
00366 //
00367 // Note : 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 }
00385 
00386 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....

Generated on Mon May 27 17:48:10 2013 for Geant4 by  doxygen 1.4.7