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

#include <G4HeatedKleinNishinaCompton.hh>

Inheritance diagram for G4HeatedKleinNishinaCompton:
G4VEmModel

Public Member Functions

 G4HeatedKleinNishinaCompton (const G4ParticleDefinition *p=0, const G4String &nam="Heated-Klein-Nishina")
 
virtual ~G4HeatedKleinNishinaCompton ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
void SetElectronTemperature (G4double t)
 
G4double GetElectronTemperature ()
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
virtual G4double ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector
< G4EmElementSelector * > * 
GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool LPMFlag () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

Protected Attributes

G4ParticleDefinitiontheGamma
 
G4ParticleDefinitiontheElectron
 
G4ParticleChangeForGammafParticleChange
 
G4double lowestGammaEnergy
 
G4double fTemperature
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Additional Inherited Members

- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Detailed Description

Definition at line 59 of file G4HeatedKleinNishinaCompton.hh.

Constructor & Destructor Documentation

G4HeatedKleinNishinaCompton::G4HeatedKleinNishinaCompton ( const G4ParticleDefinition p = 0,
const G4String nam = "Heated-Klein-Nishina" 
)
G4HeatedKleinNishinaCompton::~G4HeatedKleinNishinaCompton ( )
virtual

Definition at line 79 of file G4HeatedKleinNishinaCompton.cc.

80 {}

Member Function Documentation

G4double G4HeatedKleinNishinaCompton::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  kinEnergy,
G4double  Z,
G4double  A,
G4double  cut,
G4double  emax 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 94 of file G4HeatedKleinNishinaCompton.cc.

References test::a, test::b, python.hepunit::barn, test::c, plottest35::c1, python.hepunit::electron_mass_c2, python.hepunit::keV, and G4INCL::Math::max().

99 {
100  G4double xSection = 0.0 ;
101  if ( Z < 0.9999 ) return xSection;
102  if ( GammaEnergy < 0.01*keV ) return xSection;
103  // if ( GammaEnergy > (100.*GeV/Z) ) return xSection;
104 
105  static const G4double a = 20.0 , b = 230.0 , c = 440.0;
106 
107  static const G4double
108  d1= 2.7965e-1*barn, d2=-1.8300e-1*barn, d3= 6.7527 *barn, d4=-1.9798e+1*barn,
109  e1= 1.9756e-5*barn, e2=-1.0205e-2*barn, e3=-7.3913e-2*barn, e4= 2.7079e-2*barn,
110  f1=-3.9178e-7*barn, f2= 6.8241e-5*barn, f3= 6.0480e-5*barn, f4= 3.0274e-4*barn;
111 
112  G4double p1Z = Z*(d1 + e1*Z + f1*Z*Z), p2Z = Z*(d2 + e2*Z + f2*Z*Z),
113  p3Z = Z*(d3 + e3*Z + f3*Z*Z), p4Z = Z*(d4 + e4*Z + f4*Z*Z);
114 
115  G4double T0 = 15.0*keV;
116  if (Z < 1.5) T0 = 40.0*keV;
117 
118  G4double X = max(GammaEnergy, T0) / electron_mass_c2;
119  xSection = p1Z*std::log(1.+2.*X)/X
120  + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
121 
122  // modification for low energy. (special case for Hydrogen)
123  if (GammaEnergy < T0) {
124  G4double dT0 = 1.*keV;
125  X = (T0+dT0) / electron_mass_c2 ;
126  G4double sigma = p1Z*log(1.+2*X)/X
127  + (p2Z + p3Z*X + p4Z*X*X)/(1. + a*X + b*X*X + c*X*X*X);
128  G4double c1 = -T0*(sigma-xSection)/(xSection*dT0);
129  G4double c2 = 0.150;
130  if (Z > 1.5) c2 = 0.375-0.0556*log(Z);
131  G4double y = log(GammaEnergy/T0);
132  xSection *= exp(-y*(c1+c2*y));
133  }
134  // G4cout << "e= " << GammaEnergy << " Z= " << Z << " cross= " << xSection << G4endl;
135  return xSection;
136 }
float electron_mass_c2
Definition: hepunit.py:274
T max(const T t1, const T t2)
brief Return the largest of the two arguments
double G4double
Definition: G4Types.hh:76
tuple c1
Definition: plottest35.py:14
G4double G4HeatedKleinNishinaCompton::GetElectronTemperature ( )
inline

Definition at line 86 of file G4HeatedKleinNishinaCompton.hh.

void G4HeatedKleinNishinaCompton::Initialise ( const G4ParticleDefinition ,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 84 of file G4HeatedKleinNishinaCompton.cc.

References fParticleChange, and G4VEmModel::GetParticleChangeForGamma().

86 {
88 }
G4ParticleChangeForGamma * fParticleChange
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:121
void G4HeatedKleinNishinaCompton::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple ,
const G4DynamicParticle aDynamicGamma,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 142 of file G4HeatedKleinNishinaCompton.cc.

References CLHEP::HepLorentzVector::boost(), CLHEP::HepLorentzVector::boostVector(), DBL_MIN, python.hepunit::electron_mass_c2, fParticleChange, fStopAndKill, fTemperature, G4RandomDirection(), G4UniformRand, G4DynamicParticle::Get4Momentum(), G4VParticleChange::GetLocalEnergyDeposit(), lowestGammaEnergy, CLHEP::Hep3Vector::mag(), G4VParticleChange::ProposeLocalEnergyDeposit(), G4ParticleChangeForGamma::ProposeMomentumDirection(), G4VParticleChange::ProposeTrackStatus(), CLHEP::Hep3Vector::rotateUz(), G4ParticleChangeForGamma::SetProposedKineticEnergy(), G4INCL::DeJongSpin::shoot(), CLHEP::HepLorentzVector::t(), theElectron, python.hepunit::twopi, CLHEP::Hep3Vector::unit(), and CLHEP::HepLorentzVector::vect().

147 {
148  // The scattered gamma energy is sampled according to Klein - Nishina formula.
149  // The random number techniques of Butcher & Messel are used
150  // (Nuc Phys 20(1960),15).
151  // Note : Effects due to binding of atomic electrons are negliged.
152 
153  // We start to prepare a heated electron from Maxwell distribution.
154  // Then we try to boost to the electron rest frame and make scattering.
155  // The final step is to recover new gamma 4momentum in the lab frame
156 
157  G4double eMomentumC2 = G4RandGamma::shoot(1.5,1.);
158  eMomentumC2 *= 2*electron_mass_c2*fTemperature; // electron (pc)^2
159  G4ThreeVector eMomDir = G4RandomDirection();
160  eMomDir *= std::sqrt(eMomentumC2);
161  G4double eEnergy = std::sqrt(eMomentumC2+electron_mass_c2*electron_mass_c2);
162  G4LorentzVector electron4v = G4LorentzVector(eMomDir,eEnergy);
163  G4ThreeVector bst = electron4v.boostVector();
164 
165  G4LorentzVector gamma4v = aDynamicGamma->Get4Momentum();
166  gamma4v.boost(-bst);
167 
168  G4ThreeVector gammaMomV = gamma4v.vect();
169  G4double gamEnergy0 = gammaMomV.mag();
170 
171 
172  // G4double gamEnergy0 = aDynamicGamma->GetKineticEnergy();
173  G4double E0_m = gamEnergy0 / electron_mass_c2 ;
174 
175  // G4ThreeVector gamDirection0 = /aDynamicGamma->GetMomentumDirection();
176 
177  G4ThreeVector gamDirection0 = gammaMomV/gamEnergy0;
178 
179  // sample the energy rate of the scattered gamma in the electron rest frame
180  //
181 
182  G4double epsilon, epsilonsq, onecost, sint2, greject ;
183 
184  G4double eps0 = 1./(1. + 2.*E0_m);
185  G4double epsilon0sq = eps0*eps0;
186  G4double alpha1 = - log(eps0);
187  G4double alpha2 = 0.5*(1.- epsilon0sq);
188 
189  do
190  {
191  if ( alpha1/(alpha1+alpha2) > G4UniformRand() )
192  {
193  epsilon = exp(-alpha1*G4UniformRand()); // eps0**r
194  epsilonsq = epsilon*epsilon;
195 
196  }
197  else
198  {
199  epsilonsq = epsilon0sq + (1.- epsilon0sq)*G4UniformRand();
200  epsilon = sqrt(epsilonsq);
201  };
202 
203  onecost = (1.- epsilon)/(epsilon*E0_m);
204  sint2 = onecost*(2.-onecost);
205  greject = 1. - epsilon*sint2/(1.+ epsilonsq);
206 
207  } while (greject < G4UniformRand());
208 
209  //
210  // scattered gamma angles. ( Z - axis along the parent gamma)
211  //
212 
213  G4double cosTeta = 1. - onecost;
214  G4double sinTeta = sqrt (sint2);
215  G4double Phi = twopi * G4UniformRand();
216  G4double dirx = sinTeta*cos(Phi), diry = sinTeta*sin(Phi), dirz = cosTeta;
217 
218  //
219  // update G4VParticleChange for the scattered gamma
220  //
221 
222  G4ThreeVector gamDirection1 ( dirx,diry,dirz );
223  gamDirection1.rotateUz(gamDirection0);
224  G4double gamEnergy1 = epsilon*gamEnergy0;
225  gamDirection1 *= gamEnergy1;
226 
227  G4LorentzVector gamma4vfinal = G4LorentzVector(gamDirection1,gamEnergy1);
228 
229 
230  // kinematic of the scattered electron
231  //
232 
233  G4double eKinEnergy = gamEnergy0 - gamEnergy1;
234  G4ThreeVector eDirection = gamEnergy0*gamDirection0 - gamEnergy1*gamDirection1;
235  eDirection = eDirection.unit();
236  G4double eFinalMom = std::sqrt(eKinEnergy*(eKinEnergy+2*electron_mass_c2));
237  eDirection *= eFinalMom;
238  G4LorentzVector e4vfinal = G4LorentzVector(eDirection,gamEnergy1+electron_mass_c2);
239 
240  gamma4vfinal.boost(bst);
241  e4vfinal.boost(bst);
242 
243  gamDirection1 = gamma4vfinal.vect();
244  gamEnergy1 = gamDirection1.mag();
245  gamDirection1 /= gamEnergy1;
246 
247 
248 
249 
251 
252  if( gamEnergy1 > lowestGammaEnergy )
253  {
254  gamDirection1 /= gamEnergy1;
256  }
257  else
258  {
260  gamEnergy1 += fParticleChange->GetLocalEnergyDeposit();
262  }
263 
264  eKinEnergy = e4vfinal.t()-electron_mass_c2;
265 
266  if( eKinEnergy > DBL_MIN )
267  {
268  // create G4DynamicParticle object for the electron.
269  eDirection = e4vfinal.vect();
270  G4double eFinMomMag = eDirection.mag();
271  eDirection /= eFinMomMag;
272  G4DynamicParticle* dp = new G4DynamicParticle(theElectron,eDirection,eKinEnergy);
273  fvect->push_back(dp);
274  }
275 }
ThreeVector shoot(const G4int Ap, const G4int Af)
Hep3Vector boostVector() const
G4ParticleChangeForGamma * fParticleChange
G4ThreeVector G4RandomDirection()
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
Hep3Vector vect() const
#define G4UniformRand()
Definition: Randomize.hh:87
HepLorentzVector & boost(double, double, double)
float electron_mass_c2
Definition: hepunit.py:274
G4LorentzVector Get4Momentum() const
G4double GetLocalEnergyDeposit() const
Hep3Vector unit() const
#define DBL_MIN
Definition: templates.hh:75
void SetProposedKineticEnergy(G4double proposedKinEnergy)
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
double mag() const
CLHEP::HepLorentzVector G4LorentzVector
void G4HeatedKleinNishinaCompton::SetElectronTemperature ( G4double  t)
inline

Definition at line 85 of file G4HeatedKleinNishinaCompton.hh.

References fTemperature.

Field Documentation

G4ParticleChangeForGamma* G4HeatedKleinNishinaCompton::fParticleChange
protected
G4double G4HeatedKleinNishinaCompton::fTemperature
protected
G4double G4HeatedKleinNishinaCompton::lowestGammaEnergy
protected
G4ParticleDefinition* G4HeatedKleinNishinaCompton::theElectron
protected
G4ParticleDefinition* G4HeatedKleinNishinaCompton::theGamma
protected

Definition at line 86 of file G4HeatedKleinNishinaCompton.hh.

Referenced by G4HeatedKleinNishinaCompton().


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