Geant4-11
Public Member Functions | Protected Attributes | Private Attributes
G4DeltaAngle Class Reference

#include <G4DeltaAngle.hh>

Inheritance diagram for G4DeltaAngle:
G4VEmAngularDistribution

Public Member Functions

 G4DeltaAngle (const G4DeltaAngle &)=delete
 
 G4DeltaAngle (const G4String &name="")
 
const G4StringGetName () const
 
G4DeltaAngleoperator= (const G4DeltaAngle &right)=delete
 
virtual void PrintGeneratorInformation () const
 
G4ThreeVectorSampleDirection (const G4DynamicParticle *dp, G4double kinEnergyFinal, G4int Z, const G4Material *mat=nullptr) final
 
G4ThreeVectorSampleDirectionForShell (const G4DynamicParticle *dp, G4double kinEnergyFinal, G4int Z, G4int shellIdx, const G4Material *mat=nullptr) final
 
virtual void SamplePairDirections (const G4DynamicParticle *dp, G4double elecKinEnergy, G4double posiKinEnergy, G4ThreeVector &dirElectron, G4ThreeVector &dirPositron, G4int Z=0, const G4Material *mat=nullptr)
 
 ~G4DeltaAngle () override
 

Protected Attributes

G4ThreeVector fLocalDirection
 
G4bool fPolarisation
 

Private Attributes

const G4ParticleDefinitionfElectron
 
G4String fName
 
G4int fShellIdx
 
G4int nprob
 
std::vector< G4doubleprob
 

Detailed Description

Definition at line 55 of file G4DeltaAngle.hh.

Constructor & Destructor Documentation

◆ G4DeltaAngle() [1/2]

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

Definition at line 60 of file G4DeltaAngle.cc.

61 : G4VEmAngularDistribution("deltaVI")
62{
64 nprob = 26;
65 fShellIdx = -1;
66 prob.resize(nprob,0.0);
67}
const G4ParticleDefinition * fElectron
Definition: G4DeltaAngle.hh:79
std::vector< G4double > prob
Definition: G4DeltaAngle.hh:82
static G4Electron * Electron()
Definition: G4Electron.cc:93
G4VEmAngularDistribution(const G4String &name)

References G4Electron::Electron(), fElectron, fShellIdx, nprob, and prob.

◆ ~G4DeltaAngle()

G4DeltaAngle::~G4DeltaAngle ( )
override

Definition at line 69 of file G4DeltaAngle.cc.

70{}

◆ G4DeltaAngle() [2/2]

G4DeltaAngle::G4DeltaAngle ( const G4DeltaAngle )
delete

Member Function Documentation

◆ GetName()

const G4String & G4VEmAngularDistribution::GetName ( ) const
inlineinherited

Definition at line 111 of file G4VEmAngularDistribution.hh.

112{
113 return fName;
114}

References G4VEmAngularDistribution::fName.

◆ operator=()

G4DeltaAngle & G4DeltaAngle::operator= ( const G4DeltaAngle right)
delete

◆ PrintGeneratorInformation()

void G4VEmAngularDistribution::PrintGeneratorInformation ( ) const
virtualinherited

◆ SampleDirection()

G4ThreeVector & G4DeltaAngle::SampleDirection ( const G4DynamicParticle dp,
G4double  kinEnergyFinal,
G4int  Z,
const G4Material mat = nullptr 
)
finalvirtual

Implements G4VEmAngularDistribution.

Definition at line 82 of file G4DeltaAngle.cc.

85{
87 G4int idx = fShellIdx;
88
89 // if idx is not properly defined sample shell index
90 if(idx < 0 || idx >= nShells) {
91 if(nShells> nprob) {
92 nprob = nShells;
93 prob.resize(nprob,0.0);
94 }
95 G4double sum = 0.0;
96 for(idx=0; idx<nShells; ++idx) {
99 prob[idx] = sum;
100 }
101 sum *= G4UniformRand();
102 for(idx=0; idx<nShells; ++idx) {
103 if(sum <= prob[idx]) { break; }
104 }
105 }
107 G4double cost;
108 /*
109 G4cout << "E(keV)= " << kinEnergyFinal/keV
110 << " Ebind(keV)= " << bindingEnergy
111 << " idx= " << idx << " nShells= " << nShells << G4endl;
112 */
113 G4int n = 0;
114 G4bool isOK = false;
115 static const G4int nmax = 100;
116 do {
117 ++n;
118 // the atomic electron
120 G4double eKinEnergy = bindingEnergy*x;
121 G4double ePotEnergy = bindingEnergy*(1.0 + x);
122 G4double e = kinEnergyFinal + ePotEnergy + electron_mass_c2;
123 G4double p = sqrt((e + electron_mass_c2)*(e - electron_mass_c2));
124
125 G4double totEnergy = dp->GetTotalEnergy();
126 G4double totMomentum = dp->GetTotalMomentum();
127 if(dp->GetParticleDefinition() == fElectron) {
128 totEnergy += ePotEnergy;
129 totMomentum = sqrt((totEnergy + electron_mass_c2)
130 *(totEnergy - electron_mass_c2));
131 }
132
133 G4double eTotEnergy = eKinEnergy + electron_mass_c2;
134 G4double eTotMomentum = sqrt(eKinEnergy*(eTotEnergy + electron_mass_c2));
135 G4double costet = 2*G4UniformRand() - 1;
136 G4double sintet = sqrt((1 - costet)*(1 + costet));
137
138 cost = 1.0;
139 if(n >= nmax) {
140 /*
141 G4ExceptionDescription ed;
142 ed << "### G4DeltaAngle Warning: " << n
143 << " iterations - stop the loop with cost= 1.0 "
144 << " for " << dp->GetDefinition()->GetParticleName() << "\n"
145 << " Ekin(MeV)= " << dp->GetKineticEnergy()/MeV
146 << " Efinal(MeV)= " << kinEnergyFinal/MeV
147 << " Ebinding(MeV)= " << bindingEnergy/MeV;
148 G4Exception("G4DeltaAngle::SampleDirection","em0044",
149 JustWarning, ed,"");
150 */
151 if(0.0 == bindingEnergy) { isOK = true; }
152 bindingEnergy = 0.0;
153 }
154
155 G4double x0 = p*(totMomentum + eTotMomentum*costet);
156 /*
157 G4cout << " x0= " << x0 << " p= " << p
158 << " ptot= " << totMomentum << " pe= " << eTotMomentum
159 << " e= " << e << " totMom= " << totMomentum
160 << G4endl;
161 */
162 if(x0 > 0.0) {
163 G4double x1 = p*eTotMomentum*sintet;
164 G4double x2 = totEnergy*(eTotEnergy - e) - e*eTotEnergy
165 - totMomentum*eTotMomentum*costet + electron_mass_c2*electron_mass_c2;
166 G4double y = -x2/x0;
167 if(std::abs(y) <= 1.0) {
168 cost = -(x2 + x1*sqrt(1. - y*y))/x0;
169 if(std::abs(cost) <= 1.0) { isOK = true; }
170 else { cost = 1.0; }
171 }
172
173 /*
174 G4cout << " Ekin(MeV)= " << dp->GetKineticEnergy()
175 << " e1(keV)= " << eKinEnergy/keV
176 << " e2(keV)= " << (e - electron_mass_c2)/keV
177 << " 1-cost= " << 1-cost
178 << " x0= " << x0 << " x1= " << x1 << " x2= " << x2
179 << G4endl;
180 */
181 }
182
183 // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
184 } while(!isOK);
185
186 G4double sint = sqrt((1 - cost)*(1 + cost));
188
189 fLocalDirection.set(sint*cos(phi), sint*sin(phi), cost);
191
192 return fLocalDirection;
193}
G4double G4Log(G4double x)
Definition: G4Log.hh:226
static constexpr double twopi
Definition: G4SIunits.hh:56
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
const G4int Z[17]
#define G4UniformRand()
Definition: Randomize.hh:52
void set(double x, double y, double z)
Hep3Vector & rotateUz(const Hep3Vector &)
Definition: ThreeVector.cc:33
static G4int GetNumberOfElectrons(G4int Z, G4int SubshellNb)
static G4double GetBindingEnergy(G4int Z, G4int SubshellNb)
static G4int GetNumberOfShells(G4int Z)
const G4ThreeVector & GetMomentumDirection() const
const G4ParticleDefinition * GetParticleDefinition() const
G4double GetTotalEnergy() const
G4double GetTotalMomentum() const
G4double bindingEnergy(G4int A, G4int Z)
float electron_mass_c2
Definition: hepunit.py:273

References G4InuclSpecialFunctions::bindingEnergy(), source.hepunit::electron_mass_c2, fElectron, G4VEmAngularDistribution::fLocalDirection, fShellIdx, G4Log(), G4UniformRand, G4AtomicShells::GetBindingEnergy(), G4DynamicParticle::GetMomentumDirection(), G4AtomicShells::GetNumberOfElectrons(), G4AtomicShells::GetNumberOfShells(), G4DynamicParticle::GetParticleDefinition(), G4DynamicParticle::GetTotalEnergy(), G4DynamicParticle::GetTotalMomentum(), CLHEP::detail::n, nprob, prob, CLHEP::Hep3Vector::rotateUz(), CLHEP::Hep3Vector::set(), twopi, and Z.

Referenced by SampleDirectionForShell().

◆ SampleDirectionForShell()

G4ThreeVector & G4DeltaAngle::SampleDirectionForShell ( const G4DynamicParticle dp,
G4double  kinEnergyFinal,
G4int  Z,
G4int  shellIdx,
const G4Material mat = nullptr 
)
finalvirtual

Reimplemented from G4VEmAngularDistribution.

Definition at line 73 of file G4DeltaAngle.cc.

76{
77 fShellIdx = idx;
78 return SampleDirection(dp, kinEnergyFinal,Z, mat);
79}
G4ThreeVector & SampleDirection(const G4DynamicParticle *dp, G4double kinEnergyFinal, G4int Z, const G4Material *mat=nullptr) final
Definition: G4DeltaAngle.cc:82

References fShellIdx, SampleDirection(), and Z.

◆ SamplePairDirections()

void G4VEmAngularDistribution::SamplePairDirections ( const G4DynamicParticle dp,
G4double  elecKinEnergy,
G4double  posiKinEnergy,
G4ThreeVector dirElectron,
G4ThreeVector dirPositron,
G4int  Z = 0,
const G4Material mat = nullptr 
)
virtualinherited

Field Documentation

◆ fElectron

const G4ParticleDefinition* G4DeltaAngle::fElectron
private

Definition at line 79 of file G4DeltaAngle.hh.

Referenced by G4DeltaAngle(), and SampleDirection().

◆ fLocalDirection

G4ThreeVector G4VEmAngularDistribution::fLocalDirection
protectedinherited

◆ fName

G4String G4VEmAngularDistribution::fName
privateinherited

Definition at line 108 of file G4VEmAngularDistribution.hh.

Referenced by G4VEmAngularDistribution::GetName().

◆ fPolarisation

G4bool G4VEmAngularDistribution::fPolarisation
protectedinherited

◆ fShellIdx

G4int G4DeltaAngle::fShellIdx
private

Definition at line 81 of file G4DeltaAngle.hh.

Referenced by G4DeltaAngle(), SampleDirection(), and SampleDirectionForShell().

◆ nprob

G4int G4DeltaAngle::nprob
private

Definition at line 80 of file G4DeltaAngle.hh.

Referenced by G4DeltaAngle(), and SampleDirection().

◆ prob

std::vector<G4double> G4DeltaAngle::prob
private

Definition at line 82 of file G4DeltaAngle.hh.

Referenced by G4DeltaAngle(), and SampleDirection().


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