// // ******************************************************************** // * License and Disclaimer * // * * // * The Geant4 software is copyright of the Copyright Holders of * // * the Geant4 Collaboration. It is provided under the terms and * // * conditions of the Geant4 Software License, included in the file * // * LICENSE and available at http://cern.ch/geant4/license . These * // * include a list of copyright holders. * // * * // * Neither the authors of this software system, nor their employing * // * institutes,nor the agencies providing financial support for this * // * work make any representation or warranty, express or implied, * // * regarding this software system or assume any liability for its * // * use. Please see the license in the file LICENSE and URL above * // * for the full disclaimer and the limitation of liability. * // * * // * This code implementation is the result of the scientific and * // * technical work of the GEANT4 collaboration. * // * By using, copying, modifying or distributing the software (or * // * any work based on the software) you agree to acknowledge its * // * use in resulting scientific publications, and indicate your * // * acceptance of all terms of the Geant4 Software license. * // ******************************************************************** // // $Id: G4PolarizedAnnihilationModel.cc,v 1.8 2009/04/12 17:47:58 vnivanch Exp $ // GEANT4 tag $Name: geant4-09-03-beta-cand-01 $ // // ------------------------------------------------------------------- // // GEANT4 Class file // // // File name: G4PolarizedAnnihilationModel // // Author: Andreas Schaelicke // // Creation date: 01.05.2005 // // Modifications: // 18-07-06 use newly calculated cross sections (P. Starovoitov) // 21-08-06 update interface (A. Schaelicke) // 17-11-06 add protection agaist e+ zero energy PostStep (V.Ivanchenko) // 10-07-07 copied Initialise() method from G4eeToTwoGammaModel to provide a // local ParticleChangeForGamma object and reduce overhead // in SampleSecondaries() (A. Schaelicke) // // // Class Description: // // Implementation of polarized gamma Annihilation scattering on free electron // // ------------------------------------------------------------------- #include "G4PolarizedAnnihilationModel.hh" #include "G4PolarizationManager.hh" #include "G4PolarizationHelper.hh" #include "G4StokesVector.hh" #include "G4PolarizedAnnihilationCrossSection.hh" #include "G4ParticleChangeForGamma.hh" #include "G4TrackStatus.hh" #include "G4Gamma.hh" G4PolarizedAnnihilationModel::G4PolarizedAnnihilationModel(const G4ParticleDefinition* p, const G4String& nam) : G4eeToTwoGammaModel(p,nam),crossSectionCalculator(0),gParticleChange(0), gIsInitialised(false) { crossSectionCalculator=new G4PolarizedAnnihilationCrossSection(); } G4PolarizedAnnihilationModel::~G4PolarizedAnnihilationModel() { if (crossSectionCalculator) delete crossSectionCalculator; } //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo.... void G4PolarizedAnnihilationModel::Initialise(const G4ParticleDefinition*, const G4DataVector&) { // G4eeToTwoGammaModel::Initialise(part,dv); if(gIsInitialised) return; gParticleChange = GetParticleChangeForGamma(); gIsInitialised = true; } G4double G4PolarizedAnnihilationModel::ComputeCrossSectionPerElectron( const G4ParticleDefinition* pd, G4double kinEnergy, G4double cut, G4double emax) { G4double xs = G4eeToTwoGammaModel::ComputeCrossSectionPerElectron(pd,kinEnergy, cut,emax); G4double polzz = theBeamPolarization.z()*theTargetPolarization.z(); G4double poltt = theBeamPolarization.x()*theTargetPolarization.x() + theBeamPolarization.y()*theTargetPolarization.y(); if (polzz!=0 || poltt!=0) { G4double xval,lasym,tasym; ComputeAsymmetriesPerElectron(kinEnergy,xval,lasym,tasym); xs*=(1.+polzz*lasym+poltt*tasym); } return xs; } void G4PolarizedAnnihilationModel::ComputeAsymmetriesPerElectron(G4double ene, G4double & valueX, G4double & valueA, G4double & valueT) { // *** calculate asymmetries G4double gam = 1. + ene/electron_mass_c2; G4double xs0=crossSectionCalculator->TotalXSection(0.,1.,gam, G4StokesVector::ZERO, G4StokesVector::ZERO); G4double xsA=crossSectionCalculator->TotalXSection(0.,1.,gam, G4StokesVector::P3, G4StokesVector::P3); G4double xsT1=crossSectionCalculator->TotalXSection(0.,1.,gam, G4StokesVector::P1, G4StokesVector::P1); G4double xsT2=crossSectionCalculator->TotalXSection(0.,1.,gam, G4StokesVector::P2, G4StokesVector::P2); G4double xsT=0.5*(xsT1+xsT2); valueX=xs0; valueA=xsA/xs0-1.; valueT=xsT/xs0-1.; // G4cout<DiceEpsilon()<Initialize(epsilmax, gama, 0., theBeamPolarization, theTargetPolarization); if (crossSectionCalculator->DiceEpsilon()<0) { G4cout<<"ERROR in PolarizedAnnihilationPS::PostStepDoIt\n" <<"epsilmax DiceRoutine not appropriate ! "<DiceEpsilon()<Initialize(epsil, gama, 0., theBeamPolarization, theTargetPolarization,1); treject = crossSectionCalculator->DiceEpsilon(); treject*=epsil; if (treject>gmax || treject<0.) G4cout<<"ERROR in PolarizedAnnihilationPS::PostStepDoIt\n" <<" eps ("<trejectmax) trejectmax=treject; if (ncount>1000) { G4cout<<"WARNING in PolarizedAnnihilationPS::PostStepDoIt\n" <<"eps dicing very inefficient ="<Initialize(epsil, gama, 0., theBeamPolarization, theTargetPolarization,2); G4double gdiced =crossSectionCalculator->getVar(0); gdiced += crossSectionCalculator->getVar(3)*theBeamPolarization.p3()*theTargetPolarization.p3(); gdiced += 1.*(std::fabs(crossSectionCalculator->getVar(1)) + std::fabs(crossSectionCalculator->getVar(2)))*beamTrans*targetTrans; gdiced += 1.*std::fabs(crossSectionCalculator->getVar(4)) *(std::fabs(theBeamPolarization.p3())*targetTrans + std::fabs(theTargetPolarization.p3())*beamTrans); G4double gdist = crossSectionCalculator->getVar(0); gdist += crossSectionCalculator->getVar(3)*theBeamPolarization.p3()*theTargetPolarization.p3(); gdist += crossSectionCalculator->getVar(1)*(std::cos(phi)*theBeamPolarization.p1() + std::sin(phi)*theBeamPolarization.p2()) *(std::cos(phi)*theTargetPolarization.p1() + std::sin(phi)*theTargetPolarization.p2()); gdist += crossSectionCalculator->getVar(2)*(std::cos(phi)*theBeamPolarization.p2() - std::sin(phi)*theBeamPolarization.p1()) *(std::cos(phi)*theTargetPolarization.p2() - std::sin(phi)*theTargetPolarization.p1()); gdist += crossSectionCalculator->getVar(4) *(std::cos(phi)*theBeamPolarization.p3()*theTargetPolarization.p1() + std::cos(phi)*theBeamPolarization.p1()*theTargetPolarization.p3() + std::sin(phi)*theBeamPolarization.p3()*theTargetPolarization.p2() + std::sin(phi)*theBeamPolarization.p2()*theTargetPolarization.p3()); treject = gdist/gdiced; //G4cout<<" treject = "<Initialize(epsil,gama,phi,theBeamPolarization,theTargetPolarization,2); // ********************************************************************** Phot1Direction.rotateUz(PositDirection); // create G4DynamicParticle object for the particle1 G4DynamicParticle* aParticle1= new G4DynamicParticle (G4Gamma::Gamma(), Phot1Direction, Phot1Energy); finalGamma1Polarization=crossSectionCalculator->GetPol2(); G4double n1=finalGamma1Polarization.mag2(); if (n1>1) { G4cout<<"ERROR: PolarizedAnnihilation Polarization Vector at epsil = " <SetPolarization(finalGamma1Polarization.p1(), finalGamma1Polarization.p2(), finalGamma1Polarization.p3()); fvect->push_back(aParticle1); // ********************************************************************** G4double Eratio= Phot1Energy/Phot2Energy; G4double PositP= std::sqrt(PositKinEnergy*(PositKinEnergy+2.*electron_mass_c2)); G4ThreeVector Phot2Direction (-dirx*Eratio, -diry*Eratio, (PositP-dirz*Phot1Energy)/Phot2Energy); Phot2Direction.rotateUz(PositDirection); // create G4DynamicParticle object for the particle2 G4DynamicParticle* aParticle2= new G4DynamicParticle (G4Gamma::Gamma(), Phot2Direction, Phot2Energy); // define polarization of second final state photon finalGamma2Polarization=crossSectionCalculator->GetPol3(); G4double n2=finalGamma2Polarization.mag2(); if (n2>1) { G4cout<<"ERROR: PolarizedAnnihilation Polarization Vector at epsil = "<SetPolarization(finalGamma2Polarization.p1(), finalGamma2Polarization.p2(), finalGamma2Polarization.p3()); fvect->push_back(aParticle2); }