source: trunk/source/processes/electromagnetic/polarisation/src/G4PolarizedGammaConversionModel.cc@ 1016

Last change on this file since 1016 was 1007, checked in by garnier, 17 years ago

update to geant4.9.2

File size: 5.7 KB
Line 
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,v 1.6 2007/11/01 17:32:34 schaelic Exp $
27// GEANT4 tag $Name: geant4-09-02 $
28//
29// -------------------------------------------------------------------
30//
31// GEANT4 Class file
32//
33//
34// File name: G4PolarizedGammaConversionModel
35//
36// Author: Karim Laihem
37//
38// Creation date: 19.04.2005
39//
40// Modifications:
41// 21-08-06 Modified to fit in g4.8.1 framework (A.Schaelicke)
42// 19-03-07 Add initialisation of crossSectionCalculator (VI)
43//
44// Class Description:
45//
46// Implementation of gamma convertion to e+e- in the field of a nucleus
47// including polarization transfer
48
49// -------------------------------------------------------------------
50//
51
52//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
53//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
54
55#include "G4PolarizedGammaConversionModel.hh"
56#include "G4Electron.hh"
57#include "G4Positron.hh"
58#include "G4Gamma.hh"
59#include "Randomize.hh"
60#include "G4DataVector.hh"
61#include "G4PhysicsLogVector.hh"
62#include "G4ParticleChangeForGamma.hh"
63#include "G4PolarizedPairProductionCrossSection.hh"
64#include "G4PolarizationHelper.hh"
65
66//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
67
68
69G4PolarizedGammaConversionModel::G4PolarizedGammaConversionModel(const G4ParticleDefinition* pd,
70 const G4String& nam)
71 : G4BetheHeitlerModel(pd,nam), crossSectionCalculator(0)
72{
73}
74
75//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
76
77G4PolarizedGammaConversionModel::~G4PolarizedGammaConversionModel()
78{
79 if (crossSectionCalculator) delete crossSectionCalculator;
80}
81
82//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
83
84void G4PolarizedGammaConversionModel::Initialise(const G4ParticleDefinition* pd,
85 const G4DataVector& dv)
86{
87 G4BetheHeitlerModel::Initialise(pd,dv);
88 if (!crossSectionCalculator)
89 crossSectionCalculator = new G4PolarizedPairProductionCrossSection();
90}
91
92//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
93
94void G4PolarizedGammaConversionModel::SampleSecondaries(std::vector<G4DynamicParticle*>* vdp,
95 const G4MaterialCutsCouple* couple,
96 const G4DynamicParticle* dp,
97 G4double tmin,
98 G4double maxEnergy)
99{
100 G4BetheHeitlerModel::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
111 // determine interaction plane
112 G4ThreeVector nInteractionFrame =
113 G4PolarizationHelper::GetFrame(dp->GetMomentumDirection(),
114 (*vdp)[0]->GetMomentumDirection());
115
116 // transform polarization into interaction frame
117 beamPol.InvRotateAz(nInteractionFrame,dp->GetMomentumDirection());
118
119 // calulcate polarization transfer
120 crossSectionCalculator->SetMaterial(GetCurrentElement()->GetN(), // number of nucleons
121 GetCurrentElement()->GetZ(),
122 GetCurrentElement()->GetfCoulomb());
123 crossSectionCalculator->Initialize(gamEnergy0, lepEnergy1, sintheta,
124 beamPol, G4StokesVector::ZERO);
125
126 // deterimine final state polarization
127 G4StokesVector lep1Pol = crossSectionCalculator->GetPol2();
128 lep1Pol.RotateAz(nInteractionFrame,(*vdp)[0]->GetMomentumDirection());
129 (*vdp)[0]->SetPolarization(lep1Pol.p1(),
130 lep1Pol.p2(),
131 lep1Pol.p3());
132
133 size_t num = vdp->size();
134 if (num!=2) G4cout<<" WARNING "<<num<<" secondaries in polarized pairproduction not supported!\n";
135 for (size_t i =1; i<num; ++i) {
136 G4StokesVector lep2Pol = crossSectionCalculator->GetPol3();
137 lep2Pol.RotateAz(nInteractionFrame,(*vdp)[i]->GetMomentumDirection());
138 (*vdp)[i]->SetPolarization(lep2Pol.p1(),
139 lep2Pol.p2(),
140 lep2Pol.p3());
141
142 }
143 }
144}
145
146//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
Note: See TracBrowser for help on using the repository browser.