Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4PolarizedGammaConversionModel.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: G4PolarizedGammaConversionModel.cc 68046 2013-03-13 14:31:38Z gcosmo $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4PolarizedGammaConversionModel
34 //
35 // Author: Karim Laihem
36 //
37 // Creation date: 19.04.2005
38 //
39 // Modifications:
40 // 21-08-06 Modified to fit in g4.8.1 framework (A.Schaelicke)
41 // 19-03-07 Add initialisation of crossSectionCalculator (VI)
42 //
43 // Class Description:
44 //
45 // Implementation of gamma convertion to e+e- in the field of a nucleus
46 // including polarization transfer
47 
48 // -------------------------------------------------------------------
49 //
50 
51 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
52 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
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  : G4BetheHeitlerModel(pd,nam), crossSectionCalculator(0)
71 {
72 }
73 
74 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
75 
77 {
79 }
80 
81 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
82 
84  const G4DataVector& dv)
85 {
89 }
90 
91 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
92 
93 void G4PolarizedGammaConversionModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
94  const G4MaterialCutsCouple* couple,
95  const G4DynamicParticle* dp,
96  G4double tmin,
97  G4double maxEnergy)
98 {
99  G4BetheHeitlerModel::SampleSecondaries(vdp, couple, dp, tmin, maxEnergy);
100 
101  if(vdp && vdp->size()>0) {
102  G4double gamEnergy0 = dp->GetKineticEnergy();
103  G4double lepEnergy1 = (*vdp)[0]->GetKineticEnergy();
104  G4double sintheta = dp->GetMomentumDirection().cross((*vdp)[0]->GetMomentumDirection()).mag();
105  if (sintheta>1.) sintheta=1.;
106 
107  G4StokesVector beamPol = dp->GetPolarization();
108  beamPol.SetPhoton();
109 
110  // determine interaction plane
111  G4ThreeVector nInteractionFrame =
113  (*vdp)[0]->GetMomentumDirection());
114 
115  // transform polarization into interaction frame
116  beamPol.InvRotateAz(nInteractionFrame,dp->GetMomentumDirection());
117 
118  // calulcate polarization transfer
119  crossSectionCalculator->SetMaterial(GetCurrentElement()->GetN(), // number of nucleons
120  GetCurrentElement()->GetZ(),
121  GetCurrentElement()->GetfCoulomb());
122  crossSectionCalculator->Initialize(gamEnergy0, lepEnergy1, sintheta,
123  beamPol, G4StokesVector::ZERO);
124 
125  // deterimine final state polarization
127  lep1Pol.RotateAz(nInteractionFrame,(*vdp)[0]->GetMomentumDirection());
128  (*vdp)[0]->SetPolarization(lep1Pol.p1(),
129  lep1Pol.p2(),
130  lep1Pol.p3());
131 
132  size_t num = vdp->size();
133  if (num!=2) G4cout<<" WARNING "<<num<<" secondaries in polarized pairproduction not supported!\n";
134  for (size_t i =1; i<num; ++i) {
136  lep2Pol.RotateAz(nInteractionFrame,(*vdp)[i]->GetMomentumDirection());
137  (*vdp)[i]->SetPolarization(lep2Pol.p1(),
138  lep2Pol.p2(),
139  lep2Pol.p3());
140 
141  }
142  }
143 }
144 
145 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4double GetKineticEnergy() const
G4double p2() const
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
void SetMaterial(G4double A, G4double Z, G4double coul)
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
G4double p3() const
virtual G4StokesVector GetPol2()
G4GLOB_DLL std::ostream G4cout
G4PolarizedGammaConversionModel(const G4ParticleDefinition *p=0, const G4String &nam="polConv")
const G4ThreeVector & GetMomentumDirection() const
static G4ThreeVector GetFrame(const G4ThreeVector &, const G4ThreeVector &)
G4double p1() const
void InvRotateAz(G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
const G4ThreeVector & GetPolarization() const
virtual G4StokesVector GetPol3()
virtual void Initialise(const G4ParticleDefinition *, const G4DataVector &)
Hep3Vector cross(const Hep3Vector &) const
double G4double
Definition: G4Types.hh:76
static const G4StokesVector ZERO
double mag() const
void RotateAz(G4ThreeVector nInteractionFrame, G4ThreeVector particleDirection)
virtual void Initialize(G4double, G4double, G4double, const G4StokesVector &p0, const G4StokesVector &p1, G4int flag=0)
const G4Element * GetCurrentElement() const
Definition: G4VEmModel.hh:440