Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4PolarizedPEEffectModel.cc
Go to the documentation of this file.
1 //
2 // ********************************************************************
3 // * License and Disclaimer *
4 // * *
5 // * The Geant4 software is copyright of the Copyright Holders of *
6 // * the Geant4 Collaboration. It is provided under the terms and *
7 // * conditions of the Geant4 Software License, included in the file *
8 // * LICENSE and available at http://cern.ch/geant4/license . These *
9 // * include a list of copyright holders. *
10 // * *
11 // * Neither the authors of this software system, nor their employing *
12 // * institutes,nor the agencies providing financial support for this *
13 // * work make any representation or warranty, express or implied, *
14 // * regarding this software system or assume any liability for its *
15 // * use. Please see the license in the file LICENSE and URL above *
16 // * for the full disclaimer and the limitation of liability. *
17 // * *
18 // * This code implementation is the result of the scientific and *
19 // * technical work of the GEANT4 collaboration. *
20 // * By using, copying, modifying or distributing the software (or *
21 // * any work based on the software) you agree to acknowledge its *
22 // * use in resulting scientific publications, and indicate your *
23 // * acceptance of all terms of the Geant4 Software license. *
24 // ********************************************************************
25 //
26 // $Id: G4PolarizedPEEffectModel.cc 68046 2013-03-13 14:31:38Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4PolarizedPEEffectModel
34 //
35 // Author: Andreas Schaelicke & Karim Laihem
36 //
37 // Creation date: 22.02.2007
38 //
39 // Modifications:
40 //
41 // Class Description:
42 //
43 // Implementation of Photo electric effect
44 // including polarization transfer from circularly polarised gammas
45 
46 // -------------------------------------------------------------------
47 //
48 
49 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
50 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
51 
52 #ifndef NOIONIZATIONAS
53 
55 #include "G4Electron.hh"
56 #include "G4Positron.hh"
57 #include "G4Gamma.hh"
58 #include "Randomize.hh"
59 #include "G4DataVector.hh"
60 #include "G4PhysicsLogVector.hh"
63 #include "G4PolarizationHelper.hh"
64 
65 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
66 
67 
69  const G4String& nam)
70  : G4PEEffectFluoModel(nam),crossSectionCalculator(0),verboseLevel(0)
71 {
72 }
73 
74 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
75 
77 {
78  if (crossSectionCalculator) delete crossSectionCalculator;
79 }
80 
81 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
82 
84  const G4DataVector& dv)
85 {
87  if (!crossSectionCalculator)
88  crossSectionCalculator = new G4PolarizedPEEffectCrossSection();
89 }
90 
91 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
92 
93 void G4PolarizedPEEffectModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
94  const G4MaterialCutsCouple* couple,
95  const G4DynamicParticle* dp,
96  G4double tmin,
97  G4double maxEnergy)
98 {
99  // std::vector<G4DynamicParticle*>* vdp =
100  G4PEEffectFluoModel::SampleSecondaries(vdp,couple, dp, tmin, maxEnergy);
101 
102  if(vdp && vdp->size()>0) {
103  G4double gamEnergy0 = dp->GetKineticEnergy();
104  G4double lepEnergy1 = (*vdp)[0]->GetKineticEnergy();
105  G4double sintheta = dp->GetMomentumDirection().cross((*vdp)[0]->GetMomentumDirection()).mag();
106  if (sintheta>1.) sintheta=1.;
107 
108  G4StokesVector beamPol = dp->GetPolarization();
109  beamPol.SetPhoton();
110  // G4cout<<" beamPol "<<beamPol<<G4endl;
111 
112  // determine interaction plane
113  G4ThreeVector nInteractionFrame =
115  (*vdp)[0]->GetMomentumDirection());
116  // G4cout<<" nInteractionFrame = "<<nInteractionFrame<<G4endl;
117  if (dp->GetMomentumDirection().cross((*vdp)[0]->GetMomentumDirection()).mag()<1.e-10) {
118  // G4cout<<" nInteractionFrame not well defined "<<G4endl;
119  // G4cout<<" choosing random interaction frame "<<G4endl;
120 
122  // G4cout<<"new nInteractionFrame = "<<nInteractionFrame<<G4endl;
123  }
124  /*
125  else {
126  G4ThreeVector mom1=dp->GetMomentumDirection();
127  G4ThreeVector mom2=(*vdp)[0]->GetMomentumDirection();
128  G4cout<<" mom1 = "<<mom1<<G4endl;
129  G4cout<<" mom2 = "<<mom2<<G4endl;
130  G4ThreeVector x=mom1.cross(mom2);
131  G4cout<<" mom1 x mom2 = "<<x<<" "<<x.mag()<<G4endl;
132  G4cout<<" norm = "<<(1./x.mag()*x)<<" "<<G4endl;
133  }
134  */
135 
136  // transform polarization into interaction frame
137  beamPol.InvRotateAz(nInteractionFrame,dp->GetMomentumDirection());
138 
139  // calulcate polarization transfer
140  crossSectionCalculator->SetMaterial(GetCurrentElement()->GetN(), // number of nucleons
141  GetCurrentElement()->GetZ(),
142  GetCurrentElement()->GetfCoulomb());
143  crossSectionCalculator->Initialize(gamEnergy0, lepEnergy1, sintheta,
144  beamPol, G4StokesVector::ZERO);
145 
146  // deterimine final state polarization
147  G4StokesVector lep1Pol = crossSectionCalculator->GetPol2();
148  // G4cout<<" lepPol "<<lep1Pol<<G4endl;
149  lep1Pol.RotateAz(nInteractionFrame,(*vdp)[0]->GetMomentumDirection());
150  (*vdp)[0]->SetPolarization(lep1Pol.p1(),
151  lep1Pol.p2(),
152  lep1Pol.p3());
153 
154  // G4cout<<" lepPol "<<lep1Pol<<G4endl;
155  size_t num = vdp->size();
156  if (num!=1) G4cout<<" WARNING "<<num<<" secondaries in polarized photo electric effect not supported!\n";
157  }
158 }
159 
160 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
161 #endif //NOIONIZATIONAS
static G4ThreeVector GetRandomFrame(const G4ThreeVector &)
G4double GetKineticEnergy() const
G4double p2() const
void SetMaterial(G4double A, G4double Z, G4double coul)
G4double p3() const
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
static G4ThreeVector GetFrame(const G4ThreeVector &, const G4ThreeVector &)
G4double p1() const
G4PolarizedPEEffectModel(const G4ParticleDefinition *p=0, const G4String &nam="Polarized-PhotoElectric")
void InvRotateAz(G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
const G4ThreeVector & GetPolarization() const
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
Hep3Vector cross(const Hep3Vector &) const
virtual void Initialize(G4double aGammaE, G4double aLept0E, G4double sintheta, const G4StokesVector &beamPol, const G4StokesVector &, G4int flag=0)
double G4double
Definition: G4Types.hh:76
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
static const G4StokesVector ZERO
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
double mag() const
void RotateAz(G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
void Initialise(const G4ParticleDefinition *pd, const G4DataVector &dv)
const G4Element * GetCurrentElement() const
Definition: G4VEmModel.hh:440