Geant4  10.02.p01
MyKleinNishinaCompton.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 //
28 //
29 // $Id: MyKleinNishinaCompton.cc 86064 2014-11-07 08:49:32Z gcosmo $
30 //
31 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
32 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
33 
34 #include "MyKleinNishinaCompton.hh"
36 #include "DetectorConstruction.hh"
37 
38 #include "G4Electron.hh"
39 #include "G4Gamma.hh"
40 #include "Randomize.hh"
41 #include "G4DataVector.hh"
43 #include "G4PhysicalConstants.hh"
44 
45 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
46 
47 using namespace std;
48 
50  const G4ParticleDefinition*,
51  const G4String& nam)
52  :G4KleinNishinaCompton(0,nam), fDetector(det), fMessenger(0)
53 {
56 }
57 
58 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
59 
61 {
62  delete fMessenger;
63 }
64 
65 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
66 
68  const G4Material* mat,
69  const G4ParticleDefinition* part,
70  G4double GammaEnergy,
72 {
73  G4double xsection = G4VEmModel::CrossSectionPerVolume(mat,part,GammaEnergy);
74 
75  return xsection*fCrossSectionFactor;
76 }
77 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
78 
80  std::vector<G4DynamicParticle*>* fvect,
81  const G4MaterialCutsCouple*,
82  const G4DynamicParticle* aDynamicGamma,
83  G4double,
84  G4double)
85 {
86  // The scattered gamma energy is sampled according to Klein - Nishina formula.
87  // The random number techniques of Butcher & Messel are used
88  // (Nuc Phys 20(1960),15).
89  // Note : Effects due to binding of atomic electrons are negliged.
90 
91  G4double gamEnergy0 = aDynamicGamma->GetKineticEnergy();
92  G4double E0_m = gamEnergy0 / electron_mass_c2 ;
93 
94  G4ThreeVector gamDirection0 = aDynamicGamma->GetMomentumDirection();
95 
96  //
97  // sample the energy rate of the scattered gamma
98  //
99 
100  G4double epsilon, epsilonsq, onecost, sint2, greject ;
101 
102  G4double eps0 = 1./(1. + 2.*E0_m);
103  G4double eps0sq = eps0*eps0;
104  G4double alpha1 = - log(eps0);
105  G4double alpha2 = 0.5*(1.- eps0sq);
106 
107  do {
108  if ( alpha1/(alpha1+alpha2) > G4UniformRand() ) {
109  epsilon = exp(-alpha1*G4UniformRand()); // eps0**r
110  epsilonsq = epsilon*epsilon;
111 
112  } else {
113  epsilonsq = eps0sq + (1.- eps0sq)*G4UniformRand();
114  epsilon = sqrt(epsilonsq);
115  };
116 
117  onecost = (1.- epsilon)/(epsilon*E0_m);
118  sint2 = onecost*(2.-onecost);
119  greject = 1. - epsilon*sint2/(1.+ epsilonsq);
120 
121  } while (greject < G4UniformRand());
122 
123  //
124  // scattered gamma angles. ( Z - axis along the parent gamma)
125  //
126 
127  G4double cosTeta = 1. - onecost;
128  G4double sinTeta = sqrt (sint2);
129  G4double Phi = twopi * G4UniformRand();
130  G4double dirx = sinTeta*cos(Phi), diry = sinTeta*sin(Phi), dirz = cosTeta;
131 
132  //
133  // update G4VParticleChange for the scattered gamma
134  //
135  // beam regeneration trick : restore incident beam
136 
137  G4ThreeVector gamDirection1 ( dirx,diry,dirz );
138  gamDirection1.rotateUz(gamDirection0);
139  G4double gamEnergy1 = epsilon*gamEnergy0;
142 
143  //
144  // kinematic of the scattered electron
145  //
146 
147  G4double eKinEnergy = gamEnergy0 - gamEnergy1;
148 
149  if(eKinEnergy > DBL_MIN) {
150  G4ThreeVector eDirection
151  = gamEnergy0*gamDirection0 - gamEnergy1*gamDirection1;
152  eDirection = eDirection.unit();
153 
154  // create G4DynamicParticle object for the electron.
155  G4DynamicParticle* dp
156  = new G4DynamicParticle(theElectron,eDirection,eKinEnergy);
157  fvect->push_back(dp);
158  }
159 }
160 
161 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
162 
163 
G4ParticleChangeForGamma * fParticleChange
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:258
G4double GetKineticEnergy() const
CLHEP::Hep3Vector G4ThreeVector
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kinEnergy, G4double cut, G4double emax)
MyKleinNishinaCompton(DetectorConstruction *, const G4ParticleDefinition *p=0, const G4String &nam="myKlein-Nishina")
Definition of the MyKleinNishinaCompton class.
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
#define G4UniformRand()
Definition: Randomize.hh:97
const G4ThreeVector & GetMomentumDirection() const
static const double twopi
Definition: G4SIunits.hh:75
Definition of the MyKleinNishinaMessenger class.
#define DBL_MIN
Definition: templates.hh:75
void SetProposedKineticEnergy(G4double proposedKinEnergy)
Detector construction class to demonstrate various ways of placement.
G4ParticleDefinition * theElectron
double G4double
Definition: G4Types.hh:76
virtual void SampleSecondaries(std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
double epsilon(double density, double temperature)
MyKleinNishinaMessenger * fMessenger