G4PolarizedAnnihilationModel.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: G4PolarizedAnnihilationModel.cc 69847 2013-05-16 09:36:18Z gcosmo $
00027 //
00028 // -------------------------------------------------------------------
00029 //
00030 // GEANT4 Class file
00031 //
00032 //
00033 // File name:     G4PolarizedAnnihilationModel
00034 //
00035 // Author:        Andreas Schaelicke
00036 //
00037 // Creation date: 01.05.2005
00038 //
00039 // Modifications:
00040 // 18-07-06 use newly calculated cross sections (P. Starovoitov)
00041 // 21-08-06 update interface (A. Schaelicke)
00042 // 17-11-06 add protection agaist e+ zero energy PostStep (V.Ivanchenko)
00043 // 10-07-07 copied Initialise() method from G4eeToTwoGammaModel to provide a  
00044 //          local ParticleChangeForGamma object and reduce overhead 
00045 //          in SampleSecondaries()  (A. Schaelicke)
00046 //
00047 //
00048 // Class Description:
00049 //
00050 // Implementation of polarized gamma Annihilation scattering on free electron
00051 // 
00052 
00053 // -------------------------------------------------------------------
00054 #include "G4PolarizedAnnihilationModel.hh"
00055 #include "G4PhysicalConstants.hh"
00056 #include "G4PolarizationManager.hh"
00057 #include "G4PolarizationHelper.hh"
00058 #include "G4StokesVector.hh"
00059 #include "G4PolarizedAnnihilationCrossSection.hh"
00060 #include "G4ParticleChangeForGamma.hh"
00061 #include "G4TrackStatus.hh"
00062 #include "G4Gamma.hh"
00063 
00064 G4PolarizedAnnihilationModel::G4PolarizedAnnihilationModel(const G4ParticleDefinition* p, 
00065                                                            const G4String& nam)
00066   : G4eeToTwoGammaModel(p,nam),crossSectionCalculator(0),verboseLevel(0),gParticleChange(0),
00067     gIsInitialised(false)
00068 {
00069   crossSectionCalculator=new G4PolarizedAnnihilationCrossSection();
00070 }
00071 
00072 G4PolarizedAnnihilationModel::~G4PolarizedAnnihilationModel()
00073 {
00074   if (crossSectionCalculator) delete crossSectionCalculator;
00075 }
00076 
00077 
00078 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00079 
00080 void G4PolarizedAnnihilationModel::Initialise(const G4ParticleDefinition*,
00081                                      const G4DataVector&)
00082 {
00083   //  G4eeToTwoGammaModel::Initialise(part,dv);
00084   if(gIsInitialised) return;
00085   gParticleChange = GetParticleChangeForGamma();
00086   gIsInitialised = true;
00087 }
00088 
00089 G4double G4PolarizedAnnihilationModel::ComputeCrossSectionPerElectron(
00090                                 const G4ParticleDefinition* pd,
00091                                       G4double kinEnergy, 
00092                                       G4double cut,
00093                                       G4double emax)
00094 {
00095   G4double xs = G4eeToTwoGammaModel::ComputeCrossSectionPerElectron(pd,kinEnergy,
00096                                                                 cut,emax);
00097 
00098   G4double polzz = theBeamPolarization.z()*theTargetPolarization.z();
00099   G4double poltt = theBeamPolarization.x()*theTargetPolarization.x() 
00100                  + theBeamPolarization.y()*theTargetPolarization.y();
00101   if (polzz!=0 || poltt!=0) {
00102     G4double xval,lasym,tasym;
00103     ComputeAsymmetriesPerElectron(kinEnergy,xval,lasym,tasym);
00104     xs*=(1.+polzz*lasym+poltt*tasym);
00105   }
00106 
00107   return xs;
00108 }
00109 
00110 void G4PolarizedAnnihilationModel::ComputeAsymmetriesPerElectron(G4double ene,
00111                                                G4double & valueX,
00112                                                G4double & valueA,
00113                                                G4double & valueT)
00114 {
00115   // *** calculate asymmetries
00116   G4double gam = 1. + ene/electron_mass_c2;
00117   G4double xs0=crossSectionCalculator->TotalXSection(0.,1.,gam,
00118                                                G4StokesVector::ZERO,
00119                                                G4StokesVector::ZERO);
00120   G4double xsA=crossSectionCalculator->TotalXSection(0.,1.,gam,
00121                                                G4StokesVector::P3,
00122                                                G4StokesVector::P3);
00123   G4double xsT1=crossSectionCalculator->TotalXSection(0.,1.,gam,
00124                                                G4StokesVector::P1,
00125                                                G4StokesVector::P1);
00126   G4double xsT2=crossSectionCalculator->TotalXSection(0.,1.,gam,
00127                                                G4StokesVector::P2,
00128                                                G4StokesVector::P2);
00129   G4double xsT=0.5*(xsT1+xsT2);
00130   
00131   valueX=xs0;
00132   valueA=xsA/xs0-1.;
00133   valueT=xsT/xs0-1.;
00134   //  G4cout<<valueX<<"\t"<<valueA<<"\t"<<valueT<<"   energy = "<<gam<<G4endl;
00135   if ( (valueA < -1) || (1 < valueA)) {
00136     G4cout<< " ERROR PolarizedAnnihilationPS::ComputeAsymmetries \n";
00137     G4cout<< " something wrong in total cross section calculation (valueA)\n";
00138     G4cout<<"*********** LONG "<<valueX<<"\t"<<valueA<<"\t"<<valueT<<"   energy = "<<gam<<G4endl;
00139   }
00140   if ( (valueT < -1) || (1 < valueT)) {
00141     G4cout<< " ERROR PolarizedAnnihilationPS::ComputeAsymmetries \n";
00142     G4cout<< " something wrong in total cross section calculation (valueT)\n";
00143     G4cout<<"****** TRAN "<<valueX<<"\t"<<valueA<<"\t"<<valueT<<"   energy = "<<gam<<G4endl;
00144   }
00145 }
00146 
00147 
00148 void G4PolarizedAnnihilationModel::SampleSecondaries(std::vector<G4DynamicParticle*>* fvect,
00149                                                      const G4MaterialCutsCouple* /*couple*/,
00150                                                      const G4DynamicParticle* dp,
00151                                                      G4double /*tmin*/,
00152                                                      G4double /*maxEnergy*/) 
00153 {
00154 //   G4ParticleChangeForGamma*  gParticleChange 
00155 //     = dynamic_cast<G4ParticleChangeForGamma*>(pParticleChange);
00156   const G4Track * aTrack = gParticleChange->GetCurrentTrack();
00157 
00158   // kill primary 
00159   gParticleChange->SetProposedKineticEnergy(0.);
00160   gParticleChange->ProposeTrackStatus(fStopAndKill);
00161 
00162   // V.Ivanchenko add protection against zero kin energy
00163   G4double PositKinEnergy = dp->GetKineticEnergy();
00164 
00165   if(PositKinEnergy < DBL_MIN) {
00166 
00167     G4double cosTeta = 2.*G4UniformRand()-1.;
00168     G4double sinTeta = std::sqrt((1.0 - cosTeta)*(1.0 + cosTeta));
00169     G4double phi     = twopi * G4UniformRand();
00170     G4ThreeVector dir(sinTeta*std::cos(phi), sinTeta*std::sin(phi), cosTeta);
00171     fvect->push_back( new G4DynamicParticle(G4Gamma::Gamma(), dir, electron_mass_c2));
00172     fvect->push_back( new G4DynamicParticle(G4Gamma::Gamma(),-dir, electron_mass_c2));
00173     return;
00174   }
00175 
00176   // *** obtain and save target and beam polarization ***
00177   G4PolarizationManager * polarizationManager = G4PolarizationManager::GetInstance();
00178 
00179   // obtain polarization of the beam
00180   theBeamPolarization = aTrack->GetPolarization();
00181 
00182   // obtain polarization of the media
00183   G4VPhysicalVolume*  aPVolume  = aTrack->GetVolume();
00184   G4LogicalVolume*    aLVolume  = aPVolume->GetLogicalVolume();
00185   const G4bool targetIsPolarized = polarizationManager->IsPolarized(aLVolume);
00186   theTargetPolarization = polarizationManager->GetVolumePolarization(aLVolume);
00187 
00188   // transfer target electron polarization in frame of positron
00189   if (targetIsPolarized)
00190       theTargetPolarization.rotateUz(dp->GetMomentumDirection());
00191   
00192   G4ParticleMomentum PositDirection = dp->GetMomentumDirection();
00193 
00194   // polar asymmetry:
00195   G4double polarization = theBeamPolarization.p3()*theTargetPolarization.p3();
00196 
00197   G4double gamam1 = PositKinEnergy/electron_mass_c2;
00198   G4double gama   = gamam1+1. , gamap1 = gamam1+2.;
00199   G4double sqgrate = std::sqrt(gamam1/gamap1)/2. , sqg2m1 = std::sqrt(gamam1*gamap1);
00200 
00201   // limits of the energy sampling
00202   G4double epsilmin = 0.5 - sqgrate , epsilmax = 0.5 + sqgrate;
00203   G4double epsilqot = epsilmax/epsilmin;
00204   
00205   //
00206   // sample the energy rate of the created gammas 
00207   // note: for polarized partices, the actual dicing strategy 
00208   //       will depend on the energy, and the degree of polarization !!
00209   //
00210   G4double epsil;
00211   G4double gmax=1. + std::fabs(polarization); // crude estimate
00212 
00213   //G4bool check_range=true;
00214 
00215   crossSectionCalculator->Initialize(epsilmin, gama, 0.,  theBeamPolarization, theTargetPolarization);
00216   if (crossSectionCalculator->DiceEpsilon()<0) {
00217     G4cout<<"ERROR in PolarizedAnnihilationPS::PostStepDoIt\n"
00218           <<"epsilmin DiceRoutine not appropriate ! "<<crossSectionCalculator->DiceEpsilon()<<G4endl;
00219     //check_range=false;
00220   }
00221 
00222   crossSectionCalculator->Initialize(epsilmax, gama, 0.,  theBeamPolarization, theTargetPolarization);
00223   if (crossSectionCalculator->DiceEpsilon()<0) {
00224     G4cout<<"ERROR in PolarizedAnnihilationPS::PostStepDoIt\n"
00225           <<"epsilmax DiceRoutine not appropriate ! "<<crossSectionCalculator->DiceEpsilon()<<G4endl;
00226     //check_range=false;
00227   }
00228 
00229   G4int ncount=0;
00230   G4double trejectmax=0.;
00231   G4double treject;
00232 
00233 
00234   do {
00235     // 
00236     epsil = epsilmin*std::pow(epsilqot,G4UniformRand());
00237 
00238     crossSectionCalculator->Initialize(epsil, gama, 0., theBeamPolarization, theTargetPolarization,1);
00239 
00240     treject = crossSectionCalculator->DiceEpsilon(); 
00241     treject*=epsil;
00242 
00243     if (treject>gmax  || treject<0.) 
00244       G4cout<<"ERROR in PolarizedAnnihilationPS::PostStepDoIt\n"
00245             <<" eps ("<<epsil<<") rejection does not work properly: "<<treject<<G4endl;
00246     ++ncount;
00247     if (treject>trejectmax) trejectmax=treject;
00248     if (ncount>1000) {
00249       G4cout<<"WARNING  in PolarizedAnnihilationPS::PostStepDoIt\n"
00250             <<"eps dicing very inefficient ="<<trejectmax/gmax
00251             <<", "<<treject/gmax<<".  For secondary energy = "<<epsil<<"   "<<ncount<<G4endl;
00252       break;
00253     }
00254 
00255   } while( treject < gmax*G4UniformRand() );
00256 
00257   //
00258   // scattered Gamma angles. ( Z - axis along the parent positron)
00259   //
00260    
00261   G4double cost = (epsil*gamap1-1.)/(epsil*sqg2m1);
00262   G4double sint = std::sqrt((1.+cost)*(1.-cost));
00263   G4double phi  = 0.;
00264   G4double   beamTrans = std::sqrt(sqr(theBeamPolarization.p1()) + sqr(theBeamPolarization.p2()));
00265   G4double targetTrans = std::sqrt(sqr(theTargetPolarization.p1()) + sqr(theTargetPolarization.p2()));
00266 
00267   //  G4cout<<"phi dicing START"<<G4endl;
00268   do{
00269     phi  = twopi * G4UniformRand();
00270     crossSectionCalculator->Initialize(epsil, gama, 0., theBeamPolarization, theTargetPolarization,2);
00271 
00272     G4double gdiced =crossSectionCalculator->getVar(0);
00273     gdiced += crossSectionCalculator->getVar(3)*theBeamPolarization.p3()*theTargetPolarization.p3();
00274     gdiced += 1.*(std::fabs(crossSectionCalculator->getVar(1)) 
00275                   + std::fabs(crossSectionCalculator->getVar(2)))*beamTrans*targetTrans;
00276     gdiced += 1.*std::fabs(crossSectionCalculator->getVar(4))
00277       *(std::fabs(theBeamPolarization.p3())*targetTrans + std::fabs(theTargetPolarization.p3())*beamTrans);
00278 
00279     G4double gdist = crossSectionCalculator->getVar(0);
00280     gdist += crossSectionCalculator->getVar(3)*theBeamPolarization.p3()*theTargetPolarization.p3();
00281     gdist += crossSectionCalculator->getVar(1)*(std::cos(phi)*theBeamPolarization.p1() 
00282                                                 + std::sin(phi)*theBeamPolarization.p2())
00283                                               *(std::cos(phi)*theTargetPolarization.p1() 
00284                                                 + std::sin(phi)*theTargetPolarization.p2());
00285     gdist += crossSectionCalculator->getVar(2)*(std::cos(phi)*theBeamPolarization.p2() 
00286                                                 - std::sin(phi)*theBeamPolarization.p1())
00287                                               *(std::cos(phi)*theTargetPolarization.p2() 
00288                                                 - std::sin(phi)*theTargetPolarization.p1());
00289     gdist += crossSectionCalculator->getVar(4)
00290       *(std::cos(phi)*theBeamPolarization.p3()*theTargetPolarization.p1()
00291         + std::cos(phi)*theBeamPolarization.p1()*theTargetPolarization.p3() 
00292         + std::sin(phi)*theBeamPolarization.p3()*theTargetPolarization.p2() 
00293         + std::sin(phi)*theBeamPolarization.p2()*theTargetPolarization.p3());
00294 
00295     treject = gdist/gdiced;
00296     //G4cout<<" treject = "<<treject<<" at phi = "<<phi<<G4endl;
00297      if (treject>1.+1.e-10 || treject<0){
00298        G4cout<<"!!!ERROR in PolarizedAnnihilationPS::PostStepDoIt\n"
00299              <<" phi rejection does not work properly: "<<treject<<G4endl;
00300        G4cout<<" gdiced = "<<gdiced<<G4endl;
00301        G4cout<<" gdist = "<<gdist<<G4endl;
00302        G4cout<<" epsil = "<<epsil<<G4endl;
00303      }
00304      
00305      if (treject<1.e-3) {
00306        G4cout<<"!!!ERROR in PolarizedAnnihilationPS::PostStepDoIt\n"
00307             <<" phi rejection does not work properly: "<<treject<<"\n";
00308        G4cout<<" gdiced="<<gdiced<<"   gdist="<<gdist<<"\n";
00309        G4cout<<" epsil = "<<epsil<<G4endl;
00310      }
00311 
00312   } while( treject < G4UniformRand() );
00313   //  G4cout<<"phi dicing END"<<G4endl;
00314 
00315   G4double dirx = sint*std::cos(phi) , diry = sint*std::sin(phi) , dirz = cost;
00316 
00317   //
00318   // kinematic of the created pair
00319   //
00320   G4double TotalAvailableEnergy = PositKinEnergy + 2*electron_mass_c2;
00321   G4double Phot1Energy = epsil*TotalAvailableEnergy;
00322   G4double Phot2Energy =(1.-epsil)*TotalAvailableEnergy;
00323 
00324   // *** prepare calculation of polarization transfer ***
00325   G4ThreeVector Phot1Direction (dirx, diry, dirz);
00326 
00327   // get interaction frame
00328   G4ThreeVector  nInteractionFrame = 
00329     G4PolarizationHelper::GetFrame(PositDirection,Phot1Direction);
00330      
00331   // define proper in-plane and out-of-plane component of initial spins
00332   theBeamPolarization.InvRotateAz(nInteractionFrame,PositDirection);
00333   theTargetPolarization.InvRotateAz(nInteractionFrame,PositDirection);
00334 
00335   // calculate spin transfere matrix
00336 
00337   crossSectionCalculator->Initialize(epsil,gama,phi,theBeamPolarization,theTargetPolarization,2);
00338 
00339   // **********************************************************************
00340 
00341   Phot1Direction.rotateUz(PositDirection);   
00342   // create G4DynamicParticle object for the particle1  
00343   G4DynamicParticle* aParticle1= new G4DynamicParticle (G4Gamma::Gamma(),
00344                                                         Phot1Direction, Phot1Energy);
00345   finalGamma1Polarization=crossSectionCalculator->GetPol2();
00346   G4double n1=finalGamma1Polarization.mag2();
00347   if (n1>1) {
00348     G4cout<<"ERROR: PolarizedAnnihilation Polarization Vector at epsil = "
00349           <<epsil<<" is too large!!! \n"
00350           <<"annihi pol1= "<<finalGamma1Polarization<<", ("<<n1<<")\n";
00351     finalGamma1Polarization*=1./std::sqrt(n1);
00352   }
00353 
00354   // define polarization of first final state photon
00355   finalGamma1Polarization.SetPhoton();
00356   finalGamma1Polarization.RotateAz(nInteractionFrame,Phot1Direction);
00357   aParticle1->SetPolarization(finalGamma1Polarization.p1(),
00358                               finalGamma1Polarization.p2(),
00359                               finalGamma1Polarization.p3());
00360 
00361   fvect->push_back(aParticle1);
00362 
00363 
00364   // **********************************************************************
00365 
00366   G4double Eratio= Phot1Energy/Phot2Energy;
00367   G4double PositP= std::sqrt(PositKinEnergy*(PositKinEnergy+2.*electron_mass_c2));
00368   G4ThreeVector Phot2Direction (-dirx*Eratio, -diry*Eratio,
00369                                 (PositP-dirz*Phot1Energy)/Phot2Energy); 
00370   Phot2Direction.rotateUz(PositDirection); 
00371   // create G4DynamicParticle object for the particle2 
00372   G4DynamicParticle* aParticle2= new G4DynamicParticle (G4Gamma::Gamma(),
00373                                                         Phot2Direction, Phot2Energy);
00374 
00375   // define polarization of second final state photon
00376   finalGamma2Polarization=crossSectionCalculator->GetPol3();
00377   G4double n2=finalGamma2Polarization.mag2();
00378   if (n2>1) {
00379     G4cout<<"ERROR: PolarizedAnnihilation Polarization Vector at epsil = "<<epsil<<" is too large!!! \n";
00380     G4cout<<"annihi pol2= "<<finalGamma2Polarization<<", ("<<n2<<")\n";
00381     
00382     finalGamma2Polarization*=1./std::sqrt(n2);
00383   }
00384   finalGamma2Polarization.SetPhoton();
00385   finalGamma2Polarization.RotateAz(nInteractionFrame,Phot2Direction);
00386   aParticle2->SetPolarization(finalGamma2Polarization.p1(),
00387                               finalGamma2Polarization.p2(),
00388                               finalGamma2Polarization.p3());
00389 
00390   fvect->push_back(aParticle2);
00391 }

Generated on Mon May 27 17:49:22 2013 for Geant4 by  doxygen 1.4.7