00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025
00026
00027
00028
00029
00030
00031
00032
00033
00034
00035
00036
00037
00038
00039
00040
00041
00042
00043
00044
00045
00046
00047
00048
00049
00050
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
00066
00067
00068 G4PolarizedPEEffectModel::G4PolarizedPEEffectModel(const G4ParticleDefinition*,
00069 const G4String& nam)
00070 : G4PEEffectFluoModel(nam),crossSectionCalculator(0),verboseLevel(0)
00071 {
00072 }
00073
00074
00075
00076 G4PolarizedPEEffectModel::~G4PolarizedPEEffectModel()
00077 {
00078 if (crossSectionCalculator) delete crossSectionCalculator;
00079 }
00080
00081
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
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
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
00111
00112
00113 G4ThreeVector nInteractionFrame =
00114 G4PolarizationHelper::GetFrame(dp->GetMomentumDirection(),
00115 (*vdp)[0]->GetMomentumDirection());
00116
00117 if (dp->GetMomentumDirection().cross((*vdp)[0]->GetMomentumDirection()).mag()<1.e-10) {
00118
00119
00120
00121 nInteractionFrame = G4PolarizationHelper::GetRandomFrame(dp->GetMomentumDirection());
00122
00123 }
00124
00125
00126
00127
00128
00129
00130
00131
00132
00133
00134
00135
00136
00137 beamPol.InvRotateAz(nInteractionFrame,dp->GetMomentumDirection());
00138
00139
00140 crossSectionCalculator->SetMaterial(GetCurrentElement()->GetN(),
00141 GetCurrentElement()->GetZ(),
00142 GetCurrentElement()->GetfCoulomb());
00143 crossSectionCalculator->Initialize(gamEnergy0, lepEnergy1, sintheta,
00144 beamPol, G4StokesVector::ZERO);
00145
00146
00147 G4StokesVector lep1Pol = crossSectionCalculator->GetPol2();
00148
00149 lep1Pol.RotateAz(nInteractionFrame,(*vdp)[0]->GetMomentumDirection());
00150 (*vdp)[0]->SetPolarization(lep1Pol.p1(),
00151 lep1Pol.p2(),
00152 lep1Pol.p3());
00153
00154
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
00161 #endif //NOIONIZATIONAS