Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4DeltaAngle Class Reference

#include <G4DeltaAngle.hh>

Inheritance diagram for G4DeltaAngle:
G4VEmAngularDistribution

Public Member Functions

 G4DeltaAngle (const G4String &name="")
 
virtual ~G4DeltaAngle ()
 
virtual G4ThreeVectorSampleDirection (const G4DynamicParticle *dp, G4double kinEnergyFinal, G4int Z, const G4Material *mat=0)
 
void PrintGeneratorInformation () const
 
- Public Member Functions inherited from G4VEmAngularDistribution
 G4VEmAngularDistribution (const G4String &name)
 
virtual ~G4VEmAngularDistribution ()
 
virtual G4ThreeVectorSampleDirectionForShell (const G4DynamicParticle *dp, G4double finalTotalEnergy, G4int Z, G4int shellID, const G4Material *)
 
const G4StringGetName () const
 

Additional Inherited Members

- Protected Attributes inherited from G4VEmAngularDistribution
G4ThreeVector fLocalDirection
 

Detailed Description

Definition at line 57 of file G4DeltaAngle.hh.

Constructor & Destructor Documentation

G4DeltaAngle::G4DeltaAngle ( const G4String name = "")

Definition at line 60 of file G4DeltaAngle.cc.

References G4Electron::Electron().

61  : G4VEmAngularDistribution("deltaVI")
62 {
63  fElectron = G4Electron::Electron();
64  nprob = 26;
65  prob.resize(nprob,0.0);
66 }
G4VEmAngularDistribution(const G4String &name)
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4DeltaAngle::~G4DeltaAngle ( )
virtual

Definition at line 68 of file G4DeltaAngle.cc.

69 {}

Member Function Documentation

void G4DeltaAngle::PrintGeneratorInformation ( ) const

Definition at line 154 of file G4DeltaAngle.cc.

155 {}
G4ThreeVector & G4DeltaAngle::SampleDirection ( const G4DynamicParticle dp,
G4double  kinEnergyFinal,
G4int  Z,
const G4Material mat = 0 
)
virtual

Implements G4VEmAngularDistribution.

Definition at line 72 of file G4DeltaAngle.cc.

References CLHEP::HepLorentzVector::beta(), G4InuclSpecialFunctions::bindingEnergy(), CLHEP::HepLorentzVector::boost(), CLHEP::HepLorentzVector::boostVector(), python.hepunit::electron_mass_c2, G4VEmAngularDistribution::fLocalDirection, G4Log(), G4UniformRand, G4InuclParticleNames::gam, CLHEP::HepLorentzVector::gamma(), G4AtomicShells::GetBindingEnergy(), G4DynamicParticle::GetMomentumDirection(), G4AtomicShells::GetNumberOfElectrons(), G4AtomicShells::GetNumberOfShells(), G4DynamicParticle::GetParticleDefinition(), G4ParticleDefinition::GetPDGMass(), G4DynamicParticle::GetTotalEnergy(), G4DynamicParticle::GetTotalMomentum(), CLHEP::HepLorentzVector::mag(), CLHEP::Hep3Vector::rotateUz(), CLHEP::Hep3Vector::set(), python.hepunit::twopi, CLHEP::Hep3Vector::unit(), test::x, CLHEP::HepLorentzVector::x(), CLHEP::HepLorentzVector::y(), and CLHEP::HepLorentzVector::z().

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 }
void set(double x, double y, double z)
G4double GetTotalEnergy() const
int G4int
Definition: G4Types.hh:78
G4double GetTotalMomentum() const
#define G4UniformRand()
Definition: Randomize.hh:87
static G4double GetBindingEnergy(G4int Z, G4int SubshellNb)
const G4ThreeVector & GetMomentumDirection() const
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:72
float electron_mass_c2
Definition: hepunit.py:274
const G4ParticleDefinition * GetParticleDefinition() const
G4double G4Log(G4double x)
Definition: G4Log.hh:227
G4double GetPDGMass() const
Hep3Vector unit() const
static G4int GetNumberOfElectrons(G4int Z, G4int SubshellNb)
double G4double
Definition: G4Types.hh:76
G4double bindingEnergy(G4int A, G4int Z)
static G4int GetNumberOfShells(G4int Z)

The documentation for this class was generated from the following files: