Geant4_10
G4DeltaAngle.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: G4DeltaAngle.cc 68380 2013-03-22 18:39:29Z vnivanch $
27 //
28 // -------------------------------------------------------------------
29 //
30 // GEANT4 Class file
31 //
32 //
33 // File name: G4DeltaAngle
34 //
35 // Author: Vladimir Ivantcheko
36 //
37 // Creation date: 23 August 2013
38 //
39 // Modifications:
40 //
41 // Class Description:
42 //
43 // Delta-electron Angular Distribution Generation
44 //
45 // Class Description: End
46 //
47 // -------------------------------------------------------------------
48 //
49 
50 #include "G4DeltaAngle.hh"
51 #include "G4PhysicalConstants.hh"
52 #include "Randomize.hh"
53 #include "G4ParticleDefinition.hh"
54 #include "G4Electron.hh"
55 #include "G4AtomicShells.hh"
56 #include "G4Log.hh"
57 
58 using namespace std;
59 
61  : G4VEmAngularDistribution("deltaVI")
62 {
63  fElectron = G4Electron::Electron();
64  nprob = 26;
65  prob.resize(nprob,0.0);
66 }
67 
69 {}
70 
73  G4double kinEnergyFinal, G4int Z,
74  const G4Material*)
75 {
77  if(nShells> nprob) {
78  nprob = nShells;
79  prob.resize(nprob,0.0);
80  }
81  G4int idx;
82  G4double sum = 0.0;
83  for(idx=0; idx<nShells; ++idx) {
86  prob[idx] = sum;
87  }
88  sum *= G4UniformRand();
89  for(idx=0; idx<nShells; ++idx) {
90  if(sum <= prob[idx]) { break; }
91  }
94 
95  G4ThreeVector bst(0.0,0.0,0.0);
96  G4double cost, en, mom;
97 
98  do {
99 
100  // the atomic electron
102  G4double eKinEnergy = bindingEnergy*x;
103  G4double ePotEnergy = bindingEnergy*(1.0 + x);
104  G4double e = kinEnergyFinal + ePotEnergy + electron_mass_c2;
105 
106  G4double totEnergy = dp->GetTotalEnergy();
107  G4double totMomentum = dp->GetTotalMomentum();
108  if(dp->GetParticleDefinition() == fElectron) {
109  totEnergy += ePotEnergy;
110  totMomentum = sqrt((totEnergy + electron_mass_c2)
111  *(totEnergy - electron_mass_c2));
112  }
113 
114  G4double eTotMomentum = sqrt(eKinEnergy*(eKinEnergy + 2*electron_mass_c2));
115  G4double phi = G4UniformRand()*twopi;
116  G4double costet = 2*G4UniformRand() - 1;
117  G4double sintet = sqrt((1 - costet)*(1 + costet));
118 
119  G4LorentzVector lv0(eTotMomentum*sintet*cos(phi),
120  eTotMomentum*sintet*sin(phi),
121  eTotMomentum*costet + totMomentum,
122  eKinEnergy + electron_mass_c2 + totEnergy);
123  bst = lv0.boostVector();
124 
125  G4double m0 = lv0.mag();
126  G4double bet = lv0.beta();
127  G4double gam = lv0.gamma();
128 
129  en = 0.5*(m0*m0 - mass*mass + electron_mass_c2*electron_mass_c2)/m0;
130  mom = sqrt((en + electron_mass_c2)*(en - electron_mass_c2));
131 
132  cost= (e/gam - en)/(mom*bet);
133 
134  //G4cout << "e= " << e << " gam= " << gam << " en= " << en
135  // << " mom= " << mom << " beta= " << bet << " cost= " << cost
136  // << G4endl;
137 
138  } while(std::fabs(cost) > 1.0);
139 
140  G4double sint = sqrt((1 - cost)*(1 + cost));
141  G4double phi = twopi*G4UniformRand();
142 
143  G4LorentzVector lv1(sint*std::cos(phi)*mom, sint*std::sin(phi)*mom,
144  mom*cost, en);
145  lv1.boost(bst);
146 
147  fLocalDirection.set(lv1.x(), lv1.y(), lv1.z());
150 
151  return fLocalDirection;
152 }
153 
155 {}
void set(double x, double y, double z)
Hep3Vector boostVector() const
G4DeltaAngle(const G4String &name="")
Definition: G4DeltaAngle.cc:60
G4double GetTotalEnergy() const
tuple x
Definition: test.py:50
int G4int
Definition: G4Types.hh:78
G4double GetTotalMomentum() const
#define G4UniformRand()
Definition: Randomize.hh:87
virtual ~G4DeltaAngle()
Definition: G4DeltaAngle.cc:68
Float_t Z
Definition: plot.C:39
double mag() const
static G4double GetBindingEnergy(G4int Z, G4int SubshellNb)
const G4ThreeVector & GetMomentumDirection() const
HepLorentzVector & boost(double, double, double)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
float electron_mass_c2
Definition: hepunit.py:274
const G4ParticleDefinition * GetParticleDefinition() const
void PrintGeneratorInformation() const
G4double G4Log(G4double x)
Definition: G4Log.hh:227
G4double GetPDGMass() const
Hep3Vector unit() const
static G4int GetNumberOfElectrons(G4int Z, G4int SubshellNb)
static G4Electron * Electron()
Definition: G4Electron.cc:94
virtual G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double kinEnergyFinal, G4int Z, const G4Material *mat=0)
Definition: G4DeltaAngle.cc:72
double G4double
Definition: G4Types.hh:76
G4double bindingEnergy(G4int A, G4int Z)
static G4int GetNumberOfShells(G4int Z)