Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4HeatedKleinNishinaCompton Class Reference

#include <G4HeatedKleinNishinaCompton.hh>

Inheritance diagram for G4HeatedKleinNishinaCompton:
Collaboration diagram for G4HeatedKleinNishinaCompton:

Public Member Functions

 G4HeatedKleinNishinaCompton (const G4ParticleDefinition *p=nullptr, const G4String &nam="Heated-Klein-Nishina")
 
virtual ~G4HeatedKleinNishinaCompton ()
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) final
 
void SetElectronTemperature (G4double t)
 
G4double GetElectronTemperature ()
 
- Public Member Functions inherited from G4KleinNishinaCompton
 G4KleinNishinaCompton (const G4ParticleDefinition *p=nullptr, const G4String &nam="Klein-Nishina")
 
virtual ~G4KleinNishinaCompton ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) override
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax) override
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
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 GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerShell (const G4ParticleDefinition *, G4int Z, G4int shellIdx, G4double kinEnergy, 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 *)
 
virtual void ModelDescription (std::ostream &outFile) const
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector
< G4EmElementSelector * > * 
GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
virtual 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=nullptr)
 
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
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 

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 *)
 
- Protected Attributes inherited from G4KleinNishinaCompton
G4ParticleDefinitiontheGamma
 
G4ParticleDefinitiontheElectron
 
G4ParticleChangeForGammafParticleChange
 
G4double lowestSecondaryEnergy
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 57 of file G4HeatedKleinNishinaCompton.hh.

Constructor & Destructor Documentation

G4HeatedKleinNishinaCompton::G4HeatedKleinNishinaCompton ( const G4ParticleDefinition p = nullptr,
const G4String nam = "Heated-Klein-Nishina" 
)
explicit

Definition at line 67 of file G4HeatedKleinNishinaCompton.cc.

69  : G4KleinNishinaCompton(p, nam)
70 {
71  fTemperature = 1.0*keV;
72 }
G4KleinNishinaCompton(const G4ParticleDefinition *p=nullptr, const G4String &nam="Klein-Nishina")
static constexpr double keV
Definition: G4SIunits.hh:216
G4HeatedKleinNishinaCompton::~G4HeatedKleinNishinaCompton ( )
virtual

Definition at line 76 of file G4HeatedKleinNishinaCompton.cc.

77 {}

Member Function Documentation

G4double G4HeatedKleinNishinaCompton::GetElectronTemperature ( )
inline

Definition at line 74 of file G4HeatedKleinNishinaCompton.hh.

74 { return fTemperature; };
void G4HeatedKleinNishinaCompton::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple ,
const G4DynamicParticle aDynamicGamma,
G4double  tmin,
G4double  maxEnergy 
)
finalvirtual

Reimplemented from G4KleinNishinaCompton.

Definition at line 81 of file G4HeatedKleinNishinaCompton.cc.

86 {
87  // do nothing below the threshold
88  if(aDynamicGamma->GetKineticEnergy() <= LowEnergyLimit()) { return; }
89 
90  // The scattered gamma energy is sampled according to Klein - Nishina formula.
91  // The random number techniques of Butcher & Messel are used
92  // (Nuc Phys 20(1960),15).
93  // Note : Effects due to binding of atomic electrons are negliged.
94 
95  // We start to prepare a heated electron from Maxwell distribution.
96  // Then we try to boost to the electron rest frame and make scattering.
97  // The final step is to recover new gamma 4momentum in the lab frame
98 
99  G4double eMomentumC2 = G4RandGamma::shoot(1.5, 1.);
100  eMomentumC2 *= 2*electron_mass_c2*fTemperature; // electron (pc)^2
101  G4ThreeVector eMomDir = G4RandomDirection();
102  eMomDir *= std::sqrt(eMomentumC2);
103  G4double eEnergy = std::sqrt(eMomentumC2+electron_mass_c2*electron_mass_c2);
104  G4LorentzVector electron4v = G4LorentzVector(eMomDir,eEnergy);
105  G4ThreeVector bst = electron4v.boostVector();
106 
107  G4LorentzVector gamma4v = aDynamicGamma->Get4Momentum();
108  gamma4v.boost(-bst);
109 
110  G4ThreeVector gammaMomV = gamma4v.vect();
111  G4double gamEnergy0 = gammaMomV.mag();
112 
113 
114  // G4double gamEnergy0 = aDynamicGamma->GetKineticEnergy();
115  G4double E0_m = gamEnergy0 / electron_mass_c2 ;
116 
117  // G4ThreeVector gamDirection0 = /aDynamicGamma->GetMomentumDirection();
118 
119  G4ThreeVector gamDirection0 = gammaMomV/gamEnergy0;
120 
121  // sample the energy rate of the scattered gamma in the electron rest frame
122  //
123 
124  G4double epsilon, epsilonsq, onecost, sint2, greject ;
125 
126  G4double eps0 = 1./(1. + 2.*E0_m);
127  G4double epsilon0sq = eps0*eps0;
128  G4double alpha1 = - G4Log(eps0);
129  G4double alpha2 = 0.5*(1.- epsilon0sq);
130 
131  G4int nloop = 0;
132  do
133  {
134  ++nloop;
135  // false interaction if too many iterations
136  if(nloop > 1000) { return; }
137 
138  if ( alpha1/(alpha1+alpha2) > G4UniformRand() )
139  {
140  epsilon = G4Exp(-alpha1*G4UniformRand()); // eps0**r
141  epsilonsq = epsilon*epsilon;
142 
143  }
144  else
145  {
146  epsilonsq = epsilon0sq + (1.- epsilon0sq)*G4UniformRand();
147  epsilon = sqrt(epsilonsq);
148  };
149 
150  onecost = (1.- epsilon)/(epsilon*E0_m);
151  sint2 = onecost*(2.-onecost);
152  greject = 1. - epsilon*sint2/(1.+ epsilonsq);
153 
154  // Loop checking, 03-Aug-2015, Vladimir Ivanchenko
155  } while (greject < G4UniformRand());
156 
157  //
158  // scattered gamma angles. ( Z - axis along the parent gamma)
159  //
160 
161  G4double cosTeta = 1. - onecost;
162  G4double sinTeta = sqrt (sint2);
163  G4double Phi = twopi * G4UniformRand();
164  G4double dirx = sinTeta*cos(Phi), diry = sinTeta*sin(Phi), dirz = cosTeta;
165 
166  //
167  // update G4VParticleChange for the scattered gamma
168  //
169 
170  G4ThreeVector gamDirection1 ( dirx,diry,dirz );
171  gamDirection1.rotateUz(gamDirection0);
172  G4double gamEnergy1 = epsilon*gamEnergy0;
173  gamDirection1 *= gamEnergy1;
174 
175  G4LorentzVector gamma4vfinal = G4LorentzVector(gamDirection1,gamEnergy1);
176 
177 
178  // kinematic of the scattered electron
179  //
180 
181  G4double eKinEnergy = gamEnergy0 - gamEnergy1;
182  G4ThreeVector eDirection = gamEnergy0*gamDirection0 - gamEnergy1*gamDirection1;
183  eDirection = eDirection.unit();
184  G4double eFinalMom = std::sqrt(eKinEnergy*(eKinEnergy+2*electron_mass_c2));
185  eDirection *= eFinalMom;
186  G4LorentzVector e4vfinal = G4LorentzVector(eDirection,gamEnergy1+electron_mass_c2);
187 
188  gamma4vfinal.boost(bst);
189  e4vfinal.boost(bst);
190 
191  gamDirection1 = gamma4vfinal.vect();
192  gamEnergy1 = gamDirection1.mag();
193  gamDirection1 /= gamEnergy1;
194 
195  G4double edep = 0.0;
196  if(gamEnergy1 > lowestSecondaryEnergy) {
199  } else {
202  edep = gamEnergy1;
203  }
204 
205  //
206  // kinematic of the scattered electron
207  //
208  eKinEnergy = e4vfinal.t()-electron_mass_c2;
209 
210  if(eKinEnergy > lowestSecondaryEnergy) {
211 
212  eDirection = e4vfinal.vect().unit();
213 
214  // create G4DynamicParticle object for the electron.
215  G4DynamicParticle* dp =
216  new G4DynamicParticle(theElectron,eDirection,eKinEnergy);
217  fvect->push_back(dp);
218  } else {
219  edep += eKinEnergy;
220  }
221  // energy balance
222  if(edep > 0.0) {
224  }
225 }
ThreeVector shoot(const G4int Ap, const G4int Af)
Hep3Vector boostVector() const
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:640
G4ParticleChangeForGamma * fParticleChange
G4double GetKineticEnergy() const
G4ThreeVector G4RandomDirection()
int G4int
Definition: G4Types.hh:78
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
static constexpr double twopi
Definition: G4SIunits.hh:76
static constexpr double electron_mass_c2
Hep3Vector vect() const
#define G4UniformRand()
Definition: Randomize.hh:97
HepLorentzVector & boost(double, double, double)
G4LorentzVector Get4Momentum() const
G4double G4Log(G4double x)
Definition: G4Log.hh:230
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
Hep3Vector unit() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
G4ParticleDefinition * theElectron
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
double mag() const
double epsilon(double density, double temperature)
CLHEP::HepLorentzVector G4LorentzVector

Here is the call graph for this function:

void G4HeatedKleinNishinaCompton::SetElectronTemperature ( G4double  t)
inline

Definition at line 73 of file G4HeatedKleinNishinaCompton.hh.

73 { fTemperature = t; };

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