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

#include <G4LowEPComptonModel.hh>

+ Inheritance diagram for G4LowEPComptonModel:

Public Member Functions

 G4LowEPComptonModel (const G4ParticleDefinition *p=0, const G4String &nam="LowEPComptonModel")
 
virtual ~G4LowEPComptonModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
- 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 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 ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0., 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 *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
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)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *)
 
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
 
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 ForceBuildTable (G4bool val)
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

Protected Attributes

G4ParticleChangeForGammafParticleChange
 
- Protected Attributes inherited from G4VEmModel
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 

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 78 of file G4LowEPComptonModel.hh.

Constructor & Destructor Documentation

◆ G4LowEPComptonModel()

G4LowEPComptonModel::G4LowEPComptonModel ( const G4ParticleDefinition p = 0,
const G4String nam = "LowEPComptonModel" 
)

Definition at line 85 of file G4LowEPComptonModel.cc.

87 :G4VEmModel(nam),fParticleChange(0),isInitialised(false),
88 scatterFunctionData(0),crossSectionHandler(0),fAtomDeexcitation(0)
89{
90 lowEnergyLimit = 250 * eV;
91 highEnergyLimit = 100 * GeV;
92
93 verboseLevel=0 ;
94 // Verbosity scale:
95 // 0 = nothing
96 // 1 = warning for energy non-conservation
97 // 2 = details of energy budget
98 // 3 = calculation of cross sections, file openings, sampling of atoms
99 // 4 = entering in methods
100
101 if( verboseLevel>0 ) {
102 G4cout << "Low energy photon Compton model is constructed " << G4endl
103 << "Energy range: "
104 << lowEnergyLimit / eV << " eV - "
105 << highEnergyLimit / GeV << " GeV"
106 << G4endl;
107 }
108
109 //Mark this model as "applicable" for atomic deexcitation
111
112}
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cout
G4ParticleChangeForGamma * fParticleChange
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:641

◆ ~G4LowEPComptonModel()

G4LowEPComptonModel::~G4LowEPComptonModel ( )
virtual

Definition at line 116 of file G4LowEPComptonModel.cc.

117{
118 delete crossSectionHandler;
119 delete scatterFunctionData;
120}

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

G4double G4LowEPComptonModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  kinEnergy,
G4double  Z,
G4double  A = 0,
G4double  cut = 0,
G4double  emax = DBL_MAX 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 177 of file G4LowEPComptonModel.cc.

182{
183 if (verboseLevel > 3) {
184 G4cout << "Calling ComputeCrossSectionPerAtom() of G4LowEPComptonModel" << G4endl;
185 }
186 if (GammaEnergy < lowEnergyLimit || GammaEnergy > highEnergyLimit) { return 0.0; }
187
188 G4double cs = crossSectionHandler->FindValue(G4int(Z), GammaEnergy);
189 return cs;
190}
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
G4double FindValue(G4int Z, G4double e) const

◆ Initialise()

void G4LowEPComptonModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector cuts 
)
virtual

Implements G4VEmModel.

Definition at line 124 of file G4LowEPComptonModel.cc.

126{
127 if (verboseLevel > 2) {
128 G4cout << "Calling G4LowEPComptonModel::Initialise()" << G4endl;
129 }
130
131 if (crossSectionHandler)
132 {
133 crossSectionHandler->Clear();
134 delete crossSectionHandler;
135 }
136 delete scatterFunctionData;
137
138 // Reading of data files - all materials are read
139 crossSectionHandler = new G4CrossSectionHandler;
140 G4String crossSectionFile = "comp/ce-cs-";
141 crossSectionHandler->LoadData(crossSectionFile);
142
143 G4VDataSetAlgorithm* scatterInterpolation = new G4LogLogInterpolation;
144 G4String scatterFile = "comp/ce-sf-";
145 scatterFunctionData = new G4CompositeEMDataSet(scatterInterpolation, 1., 1.);
146 scatterFunctionData->LoadData(scatterFile);
147
148 // For Doppler broadening
149 shellData.SetOccupancyData();
150 G4String file = "/doppler/shell-doppler";
151 shellData.LoadData(file);
152
153 InitialiseElementSelectors(particle,cuts);
154
155 if (verboseLevel > 2) {
156 G4cout << "Loaded cross section files for low energy photon Compton model" << G4endl;
157 }
158
159 if(isInitialised) { return; }
160 isInitialised = true;
161
163
164 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
165
166 if( verboseLevel>0 ) {
167 G4cout << "Low energy photon Compton model is initialized " << G4endl
168 << "Energy range: "
169 << LowEnergyLimit() / eV << " eV - "
170 << HighEnergyLimit() / GeV << " GeV"
171 << G4endl;
172 }
173}
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
void SetOccupancyData()
Definition: G4ShellData.hh:70
void LoadData(const G4String &fileName)
Definition: G4ShellData.cc:234
void LoadData(const G4String &dataFile)
virtual G4bool LoadData(const G4String &fileName)=0
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:109
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:529
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:522
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:123

◆ SampleSecondaries()

void G4LowEPComptonModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle aDynamicGamma,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 199 of file G4LowEPComptonModel.cc.

203{
204
205 // The scattered gamma energy is sampled according to Klein - Nishina formula.
206 // then accepted or rejected depending on the Scattering Function multiplied
207 // by factor from Klein - Nishina formula.
208 // Expression of the angular distribution as Klein Nishina
209 // angular and energy distribution and Scattering fuctions is taken from
210 // D. E. Cullen "A simple model of photon transport" Nucl. Instr. Meth.
211 // Phys. Res. B 101 (1995). Method of sampling with form factors is different
212 // data are interpolated while in the article they are fitted.
213 // Reference to the article is from J. Stepanek New Photon, Positron
214 // and Electron Interaction Data for GEANT in Energy Range from 1 eV to 10
215 // TeV (draft).
216 // The random number techniques of Butcher & Messel are used
217 // (Nucl Phys 20(1960),15).
218
219 G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy()/MeV;
220
221 if (verboseLevel > 3) {
222 G4cout << "G4LowEPComptonModel::SampleSecondaries() E(MeV)= "
223 << photonEnergy0/MeV << " in " << couple->GetMaterial()->GetName()
224 << G4endl;
225 }
226
227 // low-energy gamma is absorpted by this process
228 if (photonEnergy0 <= lowEnergyLimit)
229 {
233 return ;
234 }
235
236 G4double e0m = photonEnergy0 / electron_mass_c2 ;
237 G4ParticleMomentum photonDirection0 = aDynamicGamma->GetMomentumDirection();
238
239 // Select randomly one element in the current material
240 const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
241 const G4Element* elm = SelectRandomAtom(couple,particle,photonEnergy0);
242 G4int Z = (G4int)elm->GetZ();
243
244 G4double LowEPCepsilon0 = 1. / (1. + 2. * e0m);
245 G4double LowEPCepsilon0Sq = LowEPCepsilon0 * LowEPCepsilon0;
246 G4double alpha1 = -std::log(LowEPCepsilon0);
247 G4double alpha2 = 0.5 * (1. - LowEPCepsilon0Sq);
248
249 G4double wlPhoton = h_Planck*c_light/photonEnergy0;
250
251 // Sample the energy of the scattered photon
252 G4double LowEPCepsilon;
253 G4double LowEPCepsilonSq;
254 G4double oneCosT;
255 G4double sinT2;
256 G4double gReject;
257
258 do
259 {
260 if ( alpha1/(alpha1+alpha2) > G4UniformRand())
261 {
262 LowEPCepsilon = std::exp(-alpha1 * G4UniformRand());
263 LowEPCepsilonSq = LowEPCepsilon * LowEPCepsilon;
264 }
265 else
266 {
267 LowEPCepsilonSq = LowEPCepsilon0Sq + (1. - LowEPCepsilon0Sq) * G4UniformRand();
268 LowEPCepsilon = std::sqrt(LowEPCepsilonSq);
269 }
270
271 oneCosT = (1. - LowEPCepsilon) / ( LowEPCepsilon * e0m);
272 sinT2 = oneCosT * (2. - oneCosT);
273 G4double x = std::sqrt(oneCosT/2.) / (wlPhoton/cm);
274 G4double scatteringFunction = scatterFunctionData->FindValue(x,Z-1);
275 gReject = (1. - LowEPCepsilon * sinT2 / (1. + LowEPCepsilonSq)) * scatteringFunction;
276
277 } while(gReject < G4UniformRand()*Z);
278
279 G4double cosTheta = 1. - oneCosT;
280 G4double sinTheta = std::sqrt(sinT2);
281 G4double phi = twopi * G4UniformRand();
282 G4double dirx = sinTheta * std::cos(phi);
283 G4double diry = sinTheta * std::sin(phi);
284 G4double dirz = cosTheta ;
285
286
287 // Scatter photon energy and Compton electron direction - Method based on:
288 // J. M. C. Brown, M. R. Dimmock, J. E. Gillam and D. M. Paganin'
289 // "A low energy bound atomic electron Compton scattering model for Geant4"
290 // TNS ISSUE, PG, 2012
291
292 // Set constants and initialize scattering parameters
293
294 G4double vel_c = 299792458;
295 G4double momentum_au_to_nat = (pi/2.0)*1.992851740*std::pow(10.,-24.);
296 G4double e_mass_kg = 9.10938188 * std::pow(10.,-31.);
297
298 G4int maxDopplerIterations = 1000;
299 G4double bindingE = 0.;
300 G4double pEIncident = photonEnergy0 ;
301 G4double pERecoil = -1.;
302 G4double eERecoil = -1.;
303 G4double e_alpha =0.;
304 G4double e_beta = 0.;
305
306 G4double CE_emission_flag = 0.;
307 G4double ePAU = -1;
308 //G4double Alpha=0;
309 G4int shellIdx = 0;
310 G4double u_temp = 0;
311 G4double cosPhiE =0;
312 G4double sinThetaE =0;
313 G4double cosThetaE =0;
314 G4int iteration = 0;
315 do{
316
317
318 // ******************************************
319 // | Determine scatter photon energy |
320 // ******************************************
321
322 do
323 {
324 iteration++;
325
326
327 // ********************************************
328 // | Sample bound electron information |
329 // ********************************************
330
331 // Select shell based on shell occupancy
332
333 shellIdx = shellData.SelectRandomShell(Z);
334 bindingE = shellData.BindingEnergy(Z,shellIdx)/MeV;
335
336 //G4cout << "New sample" << G4endl;
337
338 // Randomly sample bound electron momentum (memento: the data set is in Atomic Units)
339 ePAU = profileData.RandomSelectMomentum(Z,shellIdx);
340
341 // Convert to SI units
342
343 G4double ePSI = ePAU * momentum_au_to_nat;
344
345 //Calculate bound electron velocity and normalise to natural units
346
347 u_temp = sqrt( ((ePSI*ePSI)*(vel_c*vel_c)) / ((e_mass_kg*e_mass_kg)*(vel_c*vel_c)+(ePSI*ePSI)))/vel_c;
348
349 // Sample incident electron direction, amorphous material, to scattering photon scattering plane
350
351 e_alpha = pi*G4UniformRand();
352 e_beta = twopi*G4UniformRand();
353
354 // Total energy of system
355
356 G4double eEIncident = electron_mass_c2 / sqrt( 1 - (u_temp*u_temp));
357 G4double systemE = eEIncident + pEIncident;
358
359
360 G4double gamma_temp = 1.0 / sqrt( 1 - (u_temp*u_temp));
361 G4double numerator = gamma_temp*electron_mass_c2*(1 - u_temp * std::cos(e_alpha));
362 G4double subdenom1 = u_temp*cosTheta*std::cos(e_alpha);
363 G4double subdenom2 = u_temp*sinTheta*std::sin(e_alpha)*std::cos(e_beta);
364 G4double denominator = (1.0 - cosTheta) + (gamma_temp*electron_mass_c2*(1 - subdenom1 - subdenom2) / pEIncident);
365 pERecoil = (numerator/denominator);
366 eERecoil = systemE - pERecoil;
367 CE_emission_flag = pEIncident - pERecoil;
368 } while ( (iteration <= maxDopplerIterations) && (CE_emission_flag < bindingE));
369
370
371
372 // End of recalculation of photon energy with Doppler broadening
373
374
375
376 // *******************************************************
377 // | Determine ejected Compton electron direction |
378 // *******************************************************
379
380 // Calculate velocity of ejected Compton electron
381
382 G4double a_temp = eERecoil / electron_mass_c2;
383 G4double u_p_temp = sqrt(1 - (1 / (a_temp*a_temp)));
384
385 // Coefficients and terms from simulatenous equations
386
387 G4double sinAlpha = std::sin(e_alpha);
388 G4double cosAlpha = std::cos(e_alpha);
389 G4double sinBeta = std::sin(e_beta);
390 G4double cosBeta = std::cos(e_beta);
391
392 G4double gamma = 1.0 / sqrt(1 - (u_temp*u_temp));
393 G4double gamma_p = 1.0 / sqrt(1 - (u_p_temp*u_p_temp));
394
395 G4double var_A = pERecoil*u_p_temp*sinTheta;
396 G4double var_B = u_p_temp* (pERecoil*cosTheta-pEIncident);
397 G4double var_C = (pERecoil-pEIncident) - ( (pERecoil*pEIncident) / (gamma_p*electron_mass_c2))*(1 - cosTheta);
398
399 G4double var_D1 = gamma*electron_mass_c2*pERecoil;
400 G4double var_D2 = (1 - (u_temp*cosTheta*cosAlpha) - (u_temp*sinTheta*cosBeta*sinAlpha));
401 G4double var_D3 = ((electron_mass_c2*electron_mass_c2)*(gamma*gamma_p - 1)) - (gamma_p*electron_mass_c2*pERecoil);
402 G4double var_D = var_D1*var_D2 + var_D3;
403
404 G4double var_E1 = ((gamma*gamma_p)*(electron_mass_c2*electron_mass_c2)*(u_temp*u_p_temp)*cosAlpha);
405 G4double var_E2 = gamma_p*electron_mass_c2*pERecoil*u_p_temp*cosTheta;
406 G4double var_E = var_E1 - var_E2;
407
408 G4double var_F1 = ((gamma*gamma_p)*(electron_mass_c2*electron_mass_c2)*(u_temp*u_p_temp)*cosBeta*sinAlpha);
409 G4double var_F2 = (gamma_p*electron_mass_c2*pERecoil*u_p_temp*sinTheta);
410 G4double var_F = var_F1 - var_F2;
411
412 G4double var_G = (gamma*gamma_p)*(electron_mass_c2*electron_mass_c2)*(u_temp*u_p_temp)*sinBeta*sinAlpha;
413
414 // Two equations form a quadratic form of Wx^2 + Yx + Z = 0
415 // Coefficents and solution to quadratic
416
417 G4double var_W1 = (var_F*var_B - var_E*var_A)*(var_F*var_B - var_E*var_A);
418 G4double var_W2 = (var_G*var_G)*(var_A*var_A) + (var_G*var_G)*(var_B*var_B);
419 G4double var_W = var_W1 + var_W2;
420
421 G4double var_Y = 2.0*(((var_A*var_D-var_F*var_C)*(var_F*var_B-var_E*var_A)) - ((var_G*var_G)*var_B*var_C));
422
423 G4double var_Z1 = (var_A*var_D - var_F*var_C)*(var_A*var_D - var_F*var_C);
424 G4double var_Z2 = (var_G*var_G)*(var_C*var_C) - (var_G*var_G)*(var_A*var_A);
425 G4double var_Z = var_Z1 + var_Z2;
426
427 G4double diff = ((var_Y*var_Y)-4*var_W*var_Z);
428
429
430 // Check if diff has FPE, if so set var_Y*var_Y and 4*var_W*var_Z to six significant figures to fix
431 if (diff < 0.0)
432 {
433
434 G4double funorder = 0.0;
435 G4double sf = 6.0;
436
437 G4double diff1 = var_Y*var_Y;
438 funorder = abs(G4lrint(std::log10(static_cast<double>(abs(diff1)))+std::numeric_limits<double>::epsilon()))+sf;
439 diff1 = G4lrint(diff1*std::pow(10.0,funorder))/(1.0*std::pow(10.0,funorder));
440
441 G4double diff2 = 4*var_W*var_Z;
442 funorder = abs(G4lrint(std::log10(static_cast<double>(abs(diff2)))+std::numeric_limits<double>::epsilon()))+sf;
443 diff2 = G4lrint(diff2*std::pow(10.0,funorder))/(1.0*std::pow(10.0,funorder));
444
445 diff = diff1 -diff2;
446
447 }
448
449
450 // Plus and minus of quadratic
451
452 G4double X_p = (-var_Y + sqrt (diff))/(2*var_W);
453 G4double X_m = (-var_Y - sqrt (diff))/(2*var_W);
454
455
456 // Randomly sample one of the two possible solutions and determin theta angle of ejected Compton electron
457 G4double ThetaE = 0.;
458 G4double sol_select = G4UniformRand();
459
460 if (sol_select < 0.5)
461 {
462 ThetaE = std::acos(X_p);
463 }
464 if (sol_select > 0.5)
465 {
466 ThetaE = std::acos(X_m);
467 }
468
469 cosThetaE = std::cos(ThetaE);
470 sinThetaE = std::sin(ThetaE);
471 G4double Theta = std::acos(cosTheta);
472
473 //Calculate electron Phi
474 G4double iSinThetaE = std::sqrt(1+std::tan((pi/2.0)-ThetaE)*std::tan((pi/2.0)-ThetaE));
475 G4double iSinTheta = std::sqrt(1+std::tan((pi/2.0)-Theta)*std::tan((pi/2.0)-Theta));
476 G4double ivar_A = iSinTheta/ (pERecoil*u_p_temp);
477 // Trigs
478 cosPhiE = (var_C - var_B*cosThetaE)*(ivar_A*iSinThetaE);
479
480
481 // Check if cosPhiE has FPE, if so set to four significant figures to fix
482 if (abs(cosPhiE) > 1.0)
483 {
484 G4double funorder = 0.0;
485 G4double sf = 4.0;
486 funorder = abs(G4lrint(std::log10(static_cast<double>(abs(cosPhiE)))+std::numeric_limits<double>::epsilon()))+sf;
487 cosPhiE = G4lrint(cosPhiE*std::pow(10.0,funorder))/(1.0*std::pow(10.0,funorder));
488
489 }
490 // End of calculation of ejection Compton electron direction
491
492 //Fix for floating point errors
493
494 } while ( (iteration <= maxDopplerIterations) && (abs(cosPhiE) > 1));
495
496 // Revert to original if maximum number of iterations threshold has been reached
497
498
499 if (iteration >= maxDopplerIterations)
500 {
501 pERecoil = photonEnergy0 ;
502 bindingE = 0.;
503 dirx=0.0;
504 diry=0.0;
505 dirz=1.0;
506 }
507
508 // Set "scattered" photon direction and energy
509
510 G4ThreeVector photonDirection1(dirx,diry,dirz);
511 photonDirection1.rotateUz(photonDirection0);
512 fParticleChange->ProposeMomentumDirection(photonDirection1) ;
513
514 if (pERecoil > 0.)
515 {
517
518 // Set ejected Compton electron direction and energy
519 G4double PhiE = std::acos(cosPhiE);
520 G4double eDirX = sinThetaE * std::cos(phi+PhiE);
521 G4double eDirY = sinThetaE * std::sin(phi+PhiE);
522 G4double eDirZ = cosThetaE;
523
524 G4double eKineticEnergy = pEIncident - pERecoil - bindingE;
525
526 G4ThreeVector eDirection(eDirX,eDirY,eDirZ);
527 eDirection.rotateUz(photonDirection0);
529 eDirection,eKineticEnergy) ;
530 fvect->push_back(dp);
531
532 }
533 else
534 {
537 }
538
539
540
541
542 // sample deexcitation
543
544 if(fAtomDeexcitation && iteration < maxDopplerIterations) {
545 G4int index = couple->GetIndex();
546 if(fAtomDeexcitation->CheckDeexcitationActiveRegion(index)) {
547 size_t nbefore = fvect->size();
549 const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
550 fAtomDeexcitation->GenerateParticles(fvect, shell, Z, index);
551 size_t nafter = fvect->size();
552 if(nafter > nbefore) {
553 for (size_t i=nbefore; i<nafter; ++i) {
554 bindingE -= ((*fvect)[i])->GetKineticEnergy();
555 }
556 }
557 }
558 }
559 if(bindingE < 0.0) { bindingE = 0.0; }
561
562}
G4AtomicShellEnumerator
@ fStopAndKill
#define G4UniformRand()
Definition: Randomize.hh:53
G4double RandomSelectMomentum(G4int Z, G4int shellIndex) const
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition: G4Electron.cc:94
G4double GetZ() const
Definition: G4Element.hh:131
const G4Material * GetMaterial() const
const G4String & GetName() const
Definition: G4Material.hh:177
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4double BindingEnergy(G4int Z, G4int shellIndex) const
Definition: G4ShellData.cc:166
G4int SelectRandomShell(G4int Z) const
Definition: G4ShellData.cc:363
G4bool CheckDeexcitationActiveRegion(G4int coupleIndex)
virtual const G4AtomicShell * GetAtomicShell(G4int Z, G4AtomicShellEnumerator shell)=0
void GenerateParticles(std::vector< G4DynamicParticle * > *secVect, const G4AtomicShell *, G4int Z, G4int coupleIndex)
virtual G4double FindValue(G4double x, G4int componentId=0) const =0
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:459
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
const G4double pi
int G4lrint(double ad)
Definition: templates.hh:163

Member Data Documentation

◆ fParticleChange

G4ParticleChangeForGamma* G4LowEPComptonModel::fParticleChange
protected

Definition at line 105 of file G4LowEPComptonModel.hh.

Referenced by Initialise(), and SampleSecondaries().


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