Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4HeatedKleinNishinaCompton Class Reference

#include <G4HeatedKleinNishinaCompton.hh>

+ Inheritance 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
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)=0
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin=0.0, G4double tmax=DBL_MAX)=0
 
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 GetPartialCrossSection (const G4Material *, G4int level, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
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)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectTargetAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double logKineticEnergy, 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 *)
 
G4int SelectIsotopeNumber (const G4Element *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=nullptr)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
G4VEmModelGetTripletModel ()
 
void SetTripletModel (G4VEmModel *)
 
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) const
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetFluctuationFlag (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
void SetUseBaseMaterials (G4bool val)
 
G4bool UseBaseMaterials () 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)
 
G4VEmModeloperator= (const G4VEmModel &right)=delete
 
 G4VEmModel (const G4VEmModel &)=delete
 

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 G4MaterialpBaseMaterial
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 
G4bool lossFlucFlag
 
G4double inveplus
 
G4double pFactor
 

Detailed Description

Definition at line 56 of file G4HeatedKleinNishinaCompton.hh.

Constructor & Destructor Documentation

◆ G4HeatedKleinNishinaCompton()

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

Definition at line 66 of file G4HeatedKleinNishinaCompton.cc.

68 : G4KleinNishinaCompton(p, nam)
69{
70 fTemperature = 1.0*keV;
71}

◆ ~G4HeatedKleinNishinaCompton()

G4HeatedKleinNishinaCompton::~G4HeatedKleinNishinaCompton ( )
virtual

Definition at line 75 of file G4HeatedKleinNishinaCompton.cc.

76{}

Member Function Documentation

◆ GetElectronTemperature()

G4double G4HeatedKleinNishinaCompton::GetElectronTemperature ( )
inline

Definition at line 73 of file G4HeatedKleinNishinaCompton.hh.

73{ return fTemperature; };

◆ SampleSecondaries()

void G4HeatedKleinNishinaCompton::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple ,
const G4DynamicParticle aDynamicGamma,
G4double  tmin,
G4double  maxEnergy 
)
finalvirtual

Reimplemented from G4KleinNishinaCompton.

Definition at line 80 of file G4HeatedKleinNishinaCompton.cc.

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

◆ SetElectronTemperature()

void G4HeatedKleinNishinaCompton::SetElectronTemperature ( G4double  t)
inline

Definition at line 72 of file G4HeatedKleinNishinaCompton.hh.

72{ fTemperature = t; };

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