G4PolarizedPEEffectModel.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: G4PolarizedPEEffectModel.cc 69847 2013-05-16 09:36:18Z gcosmo $
00027 //
00028 // -------------------------------------------------------------------
00029 //
00030 // GEANT4 Class file
00031 //
00032 //
00033 // File name:     G4PolarizedPEEffectModel
00034 //
00035 // Author:        Andreas Schaelicke & Karim Laihem
00036 //
00037 // Creation date: 22.02.2007
00038 //
00039 // Modifications:
00040 //
00041 // Class Description:
00042 //
00043 // Implementation of Photo electric effect 
00044 // including polarization transfer from circularly polarised gammas
00045 
00046 // -------------------------------------------------------------------
00047 //
00048 
00049 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00050 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00051 
00052 #ifndef NOIONIZATIONAS
00053 
00054 #include "G4PolarizedPEEffectModel.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 "G4PolarizedPEEffectCrossSection.hh"
00063 #include "G4PolarizationHelper.hh"
00064 
00065 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00066 
00067 
00068 G4PolarizedPEEffectModel::G4PolarizedPEEffectModel(const G4ParticleDefinition*,
00069                                                    const G4String& nam)
00070   : G4PEEffectFluoModel(nam),crossSectionCalculator(0),verboseLevel(0)
00071 {
00072 }
00073 
00074 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00075 
00076 G4PolarizedPEEffectModel::~G4PolarizedPEEffectModel()
00077 {
00078   if (crossSectionCalculator) delete crossSectionCalculator;
00079 }
00080 
00081 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
00082 
00083 void G4PolarizedPEEffectModel::Initialise(const G4ParticleDefinition* pd,
00084                                      const G4DataVector& dv)
00085 {
00086   G4PEEffectFluoModel::Initialise(pd,dv);
00087   if (!crossSectionCalculator)
00088     crossSectionCalculator = new G4PolarizedPEEffectCrossSection();
00089 }
00090 
00091 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00092 
00093 void G4PolarizedPEEffectModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
00094                                                  const G4MaterialCutsCouple* couple,
00095                                                  const G4DynamicParticle* dp,
00096                                                  G4double tmin,
00097                                                  G4double maxEnergy)
00098 {
00099   //  std::vector<G4DynamicParticle*>* vdp = 
00100   G4PEEffectFluoModel::SampleSecondaries(vdp,couple, dp, tmin, maxEnergy);
00101 
00102   if(vdp && vdp->size()>0) {
00103     G4double gamEnergy0 = dp->GetKineticEnergy();
00104     G4double lepEnergy1 = (*vdp)[0]->GetKineticEnergy();
00105     G4double sintheta   = dp->GetMomentumDirection().cross((*vdp)[0]->GetMomentumDirection()).mag();
00106     if (sintheta>1.) sintheta=1.;
00107 
00108     G4StokesVector beamPol = dp->GetPolarization();
00109     beamPol.SetPhoton();
00110     //    G4cout<<" beamPol "<<beamPol<<G4endl;
00111 
00112     // determine interaction plane
00113     G4ThreeVector  nInteractionFrame = 
00114       G4PolarizationHelper::GetFrame(dp->GetMomentumDirection(), 
00115                                      (*vdp)[0]->GetMomentumDirection());
00116     //    G4cout<<" nInteractionFrame = "<<nInteractionFrame<<G4endl;
00117     if (dp->GetMomentumDirection().cross((*vdp)[0]->GetMomentumDirection()).mag()<1.e-10) {
00118       //      G4cout<<" nInteractionFrame not well defined "<<G4endl;
00119       //      G4cout<<" choosing random interaction frame "<<G4endl;
00120 
00121       nInteractionFrame = G4PolarizationHelper::GetRandomFrame(dp->GetMomentumDirection());
00122       //      G4cout<<"new nInteractionFrame = "<<nInteractionFrame<<G4endl;
00123     }
00124     /*
00125     else {
00126       G4ThreeVector mom1=dp->GetMomentumDirection();
00127       G4ThreeVector mom2=(*vdp)[0]->GetMomentumDirection();
00128       G4cout<<" mom1 = "<<mom1<<G4endl;
00129       G4cout<<" mom2 = "<<mom2<<G4endl;
00130       G4ThreeVector x=mom1.cross(mom2);
00131       G4cout<<" mom1 x mom2 = "<<x<<" "<<x.mag()<<G4endl;
00132       G4cout<<"  norm       = "<<(1./x.mag()*x)<<" "<<G4endl;
00133     }
00134     */
00135 
00136     // transform polarization into interaction frame
00137     beamPol.InvRotateAz(nInteractionFrame,dp->GetMomentumDirection());
00138 
00139     // calulcate polarization transfer
00140     crossSectionCalculator->SetMaterial(GetCurrentElement()->GetN(), // number of nucleons
00141                                         GetCurrentElement()->GetZ(), 
00142                                         GetCurrentElement()->GetfCoulomb());
00143     crossSectionCalculator->Initialize(gamEnergy0, lepEnergy1, sintheta,
00144                                        beamPol, G4StokesVector::ZERO);
00145 
00146     // deterimine final state polarization
00147     G4StokesVector lep1Pol = crossSectionCalculator->GetPol2();
00148     //    G4cout<<" lepPol "<<lep1Pol<<G4endl;
00149     lep1Pol.RotateAz(nInteractionFrame,(*vdp)[0]->GetMomentumDirection());
00150     (*vdp)[0]->SetPolarization(lep1Pol.p1(),
00151                                lep1Pol.p2(),
00152                                lep1Pol.p3());
00153 
00154     //    G4cout<<" lepPol "<<lep1Pol<<G4endl;
00155     size_t num = vdp->size();
00156     if (num!=1) G4cout<<" WARNING "<<num<<" secondaries in polarized photo electric effect not supported!\n"; 
00157   }
00158 }
00159 
00160 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
00161 #endif //NOIONIZATIONAS

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