G4PolarizedGammaConversionModel.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: G4PolarizedGammaConversionModel.cc 69847 2013-05-16 09:36:18Z gcosmo $
00027 //
00028 // -------------------------------------------------------------------
00029 //
00030 // GEANT4 Class file
00031 //
00032 //
00033 // File name:     G4PolarizedGammaConversionModel
00034 //
00035 // Author:        Karim Laihem
00036 //
00037 // Creation date: 19.04.2005
00038 //
00039 // Modifications:
00040 // 21-08-06 Modified to fit in g4.8.1 framework (A.Schaelicke)
00041 // 19-03-07 Add initialisation of crossSectionCalculator (VI)
00042 //
00043 // Class Description:
00044 //
00045 // Implementation of gamma convertion to e+e- in the field of a nucleus 
00046 // including polarization transfer
00047 
00048 // -------------------------------------------------------------------
00049 //
00050 
00051 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00052 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00053 
00054 #include "G4PolarizedGammaConversionModel.hh"
00055 #include "G4Electron.hh"
00056 #include "G4Positron.hh"
00057 #include "G4Gamma.hh"
00058 #include "Randomize.hh"
00059 #include "G4DataVector.hh"
00060 #include "G4PhysicsLogVector.hh"
00061 #include "G4ParticleChangeForGamma.hh"
00062 #include "G4PolarizedPairProductionCrossSection.hh"
00063 #include "G4PolarizationHelper.hh"
00064 
00065 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00066 
00067 
00068 G4PolarizedGammaConversionModel::G4PolarizedGammaConversionModel(const G4ParticleDefinition* pd,
00069                                          const G4String& nam)
00070   : G4BetheHeitlerModel(pd,nam), crossSectionCalculator(0)
00071 {
00072 }
00073 
00074 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00075 
00076 G4PolarizedGammaConversionModel::~G4PolarizedGammaConversionModel()
00077 {
00078   if (crossSectionCalculator) delete crossSectionCalculator;
00079 }
00080 
00081 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00082 
00083 void G4PolarizedGammaConversionModel::Initialise(const G4ParticleDefinition* pd,
00084                                      const G4DataVector& dv)
00085 {
00086   G4BetheHeitlerModel::Initialise(pd,dv);
00087   if (!crossSectionCalculator)
00088     crossSectionCalculator = new G4PolarizedPairProductionCrossSection();
00089 }
00090 
00091 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00092 
00093 void G4PolarizedGammaConversionModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
00094                                                         const G4MaterialCutsCouple* couple,
00095                                                         const G4DynamicParticle* dp,
00096                                                         G4double tmin,
00097                                                         G4double maxEnergy)
00098 {
00099   G4BetheHeitlerModel::SampleSecondaries(vdp, couple, dp, tmin, maxEnergy);
00100  
00101   if(vdp && vdp->size()>0) {
00102     G4double gamEnergy0 = dp->GetKineticEnergy();
00103     G4double lepEnergy1 = (*vdp)[0]->GetKineticEnergy();
00104     G4double sintheta   = dp->GetMomentumDirection().cross((*vdp)[0]->GetMomentumDirection()).mag();
00105     if (sintheta>1.) sintheta=1.;
00106 
00107     G4StokesVector beamPol = dp->GetPolarization();
00108     beamPol.SetPhoton();
00109 
00110     // determine interaction plane
00111     G4ThreeVector  nInteractionFrame = 
00112       G4PolarizationHelper::GetFrame(dp->GetMomentumDirection(), 
00113                                      (*vdp)[0]->GetMomentumDirection());
00114 
00115     // transform polarization into interaction frame
00116     beamPol.InvRotateAz(nInteractionFrame,dp->GetMomentumDirection());
00117 
00118     // calulcate polarization transfer
00119     crossSectionCalculator->SetMaterial(GetCurrentElement()->GetN(), // number of nucleons
00120                                         GetCurrentElement()->GetZ(), 
00121                                         GetCurrentElement()->GetfCoulomb());
00122     crossSectionCalculator->Initialize(gamEnergy0, lepEnergy1, sintheta,
00123                                        beamPol, G4StokesVector::ZERO);
00124 
00125     // deterimine final state polarization
00126     G4StokesVector lep1Pol = crossSectionCalculator->GetPol2();
00127     lep1Pol.RotateAz(nInteractionFrame,(*vdp)[0]->GetMomentumDirection());
00128     (*vdp)[0]->SetPolarization(lep1Pol.p1(),
00129                                lep1Pol.p2(),
00130                                lep1Pol.p3());
00131 
00132     size_t num = vdp->size();
00133     if (num!=2) G4cout<<" WARNING "<<num<<" secondaries in polarized pairproduction not supported!\n"; 
00134     for (size_t i =1; i<num; ++i) {
00135       G4StokesVector lep2Pol = crossSectionCalculator->GetPol3();
00136       lep2Pol.RotateAz(nInteractionFrame,(*vdp)[i]->GetMomentumDirection());
00137       (*vdp)[i]->SetPolarization(lep2Pol.p1(),
00138                                  lep2Pol.p2(),
00139                                  lep2Pol.p3());
00140 
00141     }
00142   }
00143 }
00144 
00145 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......

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