Geant4 11.2.2
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=nullptr, const G4String &nam="LowEPComptonModel")
 
virtual ~G4LowEPComptonModel ()
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) override
 
void InitialiseForElement (const G4ParticleDefinition *, G4int Z) override
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
G4LowEPComptonModeloperator= (const G4LowEPComptonModel &right)=delete
 
 G4LowEPComptonModel (const G4LowEPComptonModel &)=delete
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void InitialiseForMaterial (const G4ParticleDefinition *, const G4Material *)
 
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 *, const G4double &length, G4double &eloss)
 
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 FillNumberOfSecondaries (G4int &numberOfTriplets, G4int &numberOfRecoil)
 
virtual void ModelDescription (std::ostream &outFile) const
 
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)
 
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)
 
const G4ElementGetCurrentElement (const G4Material *mat=nullptr) const
 
G4int SelectRandomAtomNumber (const G4Material *) const
 
const G4IsotopeGetCurrentIsotope (const G4Element *elm=nullptr) const
 
G4int SelectIsotopeNumber (const G4Element *) const
 
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 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 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 *)
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 
void SetLPMFlag (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 G4VEmModel
G4ElementDatafElementData = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4PhysicsTablexSectionTable = nullptr
 
const G4MaterialpBaseMaterial = nullptr
 
const std::vector< G4double > * theDensityFactor = nullptr
 
const std::vector< G4int > * theDensityIdx = nullptr
 
G4double inveplus
 
G4double pFactor = 1.0
 
std::size_t currentCoupleIndex = 0
 
std::size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

Detailed Description

Definition at line 81 of file G4LowEPComptonModel.hh.

Constructor & Destructor Documentation

◆ G4LowEPComptonModel() [1/2]

G4LowEPComptonModel::G4LowEPComptonModel ( const G4ParticleDefinition * p = nullptr,
const G4String & nam = "LowEPComptonModel" )
explicit

Definition at line 97 of file G4LowEPComptonModel.cc.

99 : G4VEmModel(nam),isInitialised(false)
100{
101 verboseLevel=1 ;
102 // Verbosity scale:
103 // 0 = nothing
104 // 1 = warning for energy non-conservation
105 // 2 = details of energy budget
106 // 3 = calculation of cross sections, file openings, sampling of atoms
107 // 4 = entering in methods
108
109 if( verboseLevel>1 ) {
110 G4cout << "Low energy photon Compton model is constructed " << G4endl;
111 }
112
113 //Mark this model as "applicable" for atomic deexcitation
115
116 fParticleChange = nullptr;
117 fAtomDeexcitation = nullptr;
118}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void SetDeexcitationFlag(G4bool val)
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67

◆ ~G4LowEPComptonModel()

G4LowEPComptonModel::~G4LowEPComptonModel ( )
virtual

Definition at line 122 of file G4LowEPComptonModel.cc.

123{
124 if(IsMaster()) {
125 delete shellData;
126 shellData = nullptr;
127 delete profileData;
128 profileData = nullptr;
129 }
130}
G4bool IsMaster() const

◆ G4LowEPComptonModel() [2/2]

G4LowEPComptonModel::G4LowEPComptonModel ( const G4LowEPComptonModel & )
delete

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 258 of file G4LowEPComptonModel.cc.

262{
263 if (verboseLevel > 3) {
264 G4cout << "G4LowEPComptonModel::ComputeCrossSectionPerAtom()"
265 << G4endl;
266 }
267 G4double cs = 0.0;
268
269 if (GammaEnergy < LowEnergyLimit()) { return 0.0; }
270
271 G4int intZ = G4lrint(Z);
272 if(intZ < 1 || intZ > maxZ) { return cs; }
273
274 G4PhysicsFreeVector* pv = data[intZ];
275
276 // if element was not initialised
277 // do initialisation safely for MT mode
278 if(!pv)
279 {
280 InitialiseForElement(0, intZ);
281 pv = data[intZ];
282 if(!pv) { return cs; }
283 }
284
285 G4int n = G4int(pv->GetVectorLength() - 1);
286 G4double e1 = pv->Energy(0);
287 G4double e2 = pv->Energy(n);
288
289 if(GammaEnergy <= e1) { cs = GammaEnergy/(e1*e1)*pv->Value(e1); }
290 else if(GammaEnergy <= e2) { cs = pv->Value(GammaEnergy)/GammaEnergy; }
291 else if(GammaEnergy > e2) { cs = pv->Value(e2)/GammaEnergy; }
292
293 return cs;
294}
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
void InitialiseForElement(const G4ParticleDefinition *, G4int Z) override
G4double Energy(const std::size_t index) const
G4double Value(const G4double energy, std::size_t &lastidx) const
std::size_t GetVectorLength() const
G4double LowEnergyLimit() const
int G4lrint(double ad)
Definition templates.hh:134

◆ Initialise()

void G4LowEPComptonModel::Initialise ( const G4ParticleDefinition * particle,
const G4DataVector & cuts )
overridevirtual

Implements G4VEmModel.

Definition at line 134 of file G4LowEPComptonModel.cc.

136{
137 if (verboseLevel > 1) {
138 G4cout << "Calling G4LowEPComptonModel::Initialise()" << G4endl;
139 }
140
141 // Initialise element selector
142 if(IsMaster()) {
143
144 // Access to elements
145 const char* path = G4FindDataDir("G4LEDATA");
146
147 G4ProductionCutsTable* theCoupleTable =
149 G4int numOfCouples = (G4int)theCoupleTable->GetTableSize();
150
151 for(G4int i=0; i<numOfCouples; ++i) {
152 const G4Material* material =
153 theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
154 const G4ElementVector* theElementVector = material->GetElementVector();
155 std::size_t nelm = material->GetNumberOfElements();
156
157 for (std::size_t j=0; j<nelm; ++j) {
158 G4int Z = G4lrint((*theElementVector)[j]->GetZ());
159 if(Z < 1) { Z = 1; }
160 else if(Z > maxZ){ Z = maxZ; }
161
162 if( (!data[Z]) ) { ReadData(Z, path); }
163 }
164 }
165
166 // For Doppler broadening
167 if(!shellData) {
168 shellData = new G4ShellData();
169 shellData->SetOccupancyData();
170 G4String file = "/doppler/shell-doppler";
171 shellData->LoadData(file);
172 }
173 if(!profileData) { profileData = new G4DopplerProfile(); }
174
175 InitialiseElementSelectors(particle, cuts);
176 }
177
178 if (verboseLevel > 2) {
179 G4cout << "Loaded cross section files" << G4endl;
180 }
181
182 if( verboseLevel>1 ) {
183 G4cout << "G4LowEPComptonModel is initialized " << G4endl
184 << "Energy range: "
185 << LowEnergyLimit() / eV << " eV - "
186 << HighEnergyLimit() / GeV << " GeV"
187 << G4endl;
188 }
189
190 if(isInitialised) { return; }
191
192 fParticleChange = GetParticleChangeForGamma();
193 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
194 isInitialised = true;
195}
std::vector< const G4Element * > G4ElementVector
const char * G4FindDataDir(const char *)
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
const G4Material * GetMaterial() const
const G4ElementVector * GetElementVector() const
std::size_t GetNumberOfElements() const
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
void SetOccupancyData()
void LoadData(const G4String &fileName)
G4ParticleChangeForGamma * GetParticleChangeForGamma()
G4double HighEnergyLimit() const
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)

◆ InitialiseForElement()

void G4LowEPComptonModel::InitialiseForElement ( const G4ParticleDefinition * ,
G4int Z )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 670 of file G4LowEPComptonModel.cc.

672{
673 G4AutoLock l(&LowEPComptonModelMutex);
674 if(!data[Z]) { ReadData(Z); }
675 l.unlock();
676}

Referenced by ComputeCrossSectionPerAtom().

◆ InitialiseLocal()

void G4LowEPComptonModel::InitialiseLocal ( const G4ParticleDefinition * ,
G4VEmModel * masterModel )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 199 of file G4LowEPComptonModel.cc.

201{
203}
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
std::vector< G4EmElementSelector * > * GetElementSelectors()

◆ operator=()

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 298 of file G4LowEPComptonModel.cc.

302{
303
304 // The scattered gamma energy is sampled according to Klein - Nishina formula.
305 // then accepted or rejected depending on the Scattering Function multiplied
306 // by factor from Klein - Nishina formula.
307 // Expression of the angular distribution as Klein Nishina
308 // angular and energy distribution and Scattering fuctions is taken from
309 // D. E. Cullen "A simple model of photon transport" Nucl. Instr. Meth.
310 // Phys. Res. B 101 (1995). Method of sampling with form factors is different
311 // data are interpolated while in the article they are fitted.
312 // Reference to the article is from J. Stepanek New Photon, Positron
313 // and Electron Interaction Data for GEANT in Energy Range from 1 eV to 10
314 // TeV (draft).
315 // The random number techniques of Butcher & Messel are used
316 // (Nucl Phys 20(1960),15).
317
318 G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy()/MeV;
319
320 if (verboseLevel > 3) {
321 G4cout << "G4LowEPComptonModel::SampleSecondaries() E(MeV)= "
322 << photonEnergy0/MeV << " in " << couple->GetMaterial()->GetName()
323 << G4endl;
324 }
325 // do nothing below the threshold
326 // should never get here because the XS is zero below the limit
327 if (photonEnergy0 < LowEnergyLimit())
328 return ;
329
330 G4double e0m = photonEnergy0 / electron_mass_c2 ;
331 G4ParticleMomentum photonDirection0 = aDynamicGamma->GetMomentumDirection();
332
333 // Select randomly one element in the current material
334 const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
335 const G4Element* elm = SelectRandomAtom(couple,particle,photonEnergy0);
336 G4int Z = (G4int)elm->GetZ();
337
338 G4double LowEPCepsilon0 = 1. / (1. + 2. * e0m);
339 G4double LowEPCepsilon0Sq = LowEPCepsilon0 * LowEPCepsilon0;
340 G4double alpha1 = -std::log(LowEPCepsilon0);
341 G4double alpha2 = 0.5 * (1. - LowEPCepsilon0Sq);
342
343 G4double wlPhoton = h_Planck*c_light/photonEnergy0;
344
345 // Sample the energy of the scattered photon
346 G4double LowEPCepsilon;
347 G4double LowEPCepsilonSq;
348 G4double oneCosT;
349 G4double sinT2;
350 G4double gReject;
351
352 if (verboseLevel > 3) {
353 G4cout << "Started loop to sample gamma energy" << G4endl;
354 }
355
356 do
357 {
358 if ( alpha1/(alpha1+alpha2) > G4UniformRand())
359 {
360 LowEPCepsilon = G4Exp(-alpha1 * G4UniformRand());
361 LowEPCepsilonSq = LowEPCepsilon * LowEPCepsilon;
362 }
363 else
364 {
365 LowEPCepsilonSq = LowEPCepsilon0Sq + (1. - LowEPCepsilon0Sq) * G4UniformRand();
366 LowEPCepsilon = std::sqrt(LowEPCepsilonSq);
367 }
368
369 oneCosT = (1. - LowEPCepsilon) / ( LowEPCepsilon * e0m);
370 sinT2 = oneCosT * (2. - oneCosT);
371 G4double x = std::sqrt(oneCosT/2.) / (wlPhoton/cm);
372 G4double scatteringFunction = ComputeScatteringFunction(x, Z);
373 gReject = (1. - LowEPCepsilon * sinT2 / (1. + LowEPCepsilonSq)) * scatteringFunction;
374
375 } while(gReject < G4UniformRand()*Z);
376
377 G4double cosTheta = 1. - oneCosT;
378 G4double sinTheta = std::sqrt(sinT2);
379 G4double phi = twopi * G4UniformRand();
380 G4double dirx = sinTheta * std::cos(phi);
381 G4double diry = sinTheta * std::sin(phi);
382 G4double dirz = cosTheta ;
383
384 // Scatter photon energy and Compton electron direction - Method based on:
385 // J. M. C. Brown, M. R. Dimmock, J. E. Gillam and D. M. Paganin'
386 // "A low energy bound atomic electron Compton scattering model for Geant4"
387 // NIMB, Vol. 338, 77-88, 2014.
388
389 // Set constants and initialize scattering parameters
390 const G4double vel_c = c_light / (m/s);
391 const G4double momentum_au_to_nat = halfpi* hbar_Planck / Bohr_radius / (kg*m/s);
392 const G4double e_mass_kg = electron_mass_c2 / c_squared / kg ;
393
394 const G4int maxDopplerIterations = 1000;
395 G4double bindingE = 0.;
396 G4double pEIncident = photonEnergy0 ;
397 G4double pERecoil = -1.;
398 G4double eERecoil = -1.;
399 G4double e_alpha =0.;
400 G4double e_beta = 0.;
401
402 G4double CE_emission_flag = 0.;
403 G4double ePAU = -1;
404 G4int shellIdx = 0;
405 G4double u_temp = 0;
406 G4double cosPhiE =0;
407 G4double sinThetaE =0;
408 G4double cosThetaE =0;
409 G4int iteration = 0;
410
411 if (verboseLevel > 3) {
412 G4cout << "Started loop to sample photon energy and electron direction" << G4endl;
413 }
414
415 do{
416 // ******************************************
417 // | Determine scatter photon energy |
418 // ******************************************
419 do
420 {
421 iteration++;
422 // ********************************************
423 // | Sample bound electron information |
424 // ********************************************
425
426 // Select shell based on shell occupancy
427 shellIdx = shellData->SelectRandomShell(Z);
428 bindingE = shellData->BindingEnergy(Z,shellIdx)/MeV;
429
430 // Randomly sample bound electron momentum (memento: the data set is in Atomic Units)
431 ePAU = profileData->RandomSelectMomentum(Z,shellIdx);
432
433 // Convert to SI units
434 G4double ePSI = ePAU * momentum_au_to_nat;
435
436 //Calculate bound electron velocity and normalise to natural units
437 u_temp = sqrt( ((ePSI*ePSI)*(vel_c*vel_c)) / ((e_mass_kg*e_mass_kg)*(vel_c*vel_c)+(ePSI*ePSI)) )/vel_c;
438
439 // Sample incident electron direction, amorphous material, to scattering photon scattering plane
440 e_alpha = pi*G4UniformRand();
441 e_beta = twopi*G4UniformRand();
442
443 // Total energy of system
444
445 G4double eEIncident = electron_mass_c2 / sqrt( 1 - (u_temp*u_temp));
446 G4double systemE = eEIncident + pEIncident;
447 G4double gamma_temp = 1.0 / sqrt( 1 - (u_temp*u_temp));
448 G4double numerator = gamma_temp*electron_mass_c2*(1 - u_temp * std::cos(e_alpha));
449 G4double subdenom1 = u_temp*cosTheta*std::cos(e_alpha);
450 G4double subdenom2 = u_temp*sinTheta*std::sin(e_alpha)*std::cos(e_beta);
451 G4double denominator = (1.0 - cosTheta) + (gamma_temp*electron_mass_c2*(1 - subdenom1 - subdenom2) / pEIncident);
452 pERecoil = (numerator/denominator);
453 eERecoil = systemE - pERecoil;
454 CE_emission_flag = pEIncident - pERecoil;
455 } while ( (iteration <= maxDopplerIterations) && (CE_emission_flag < bindingE));
456
457 // End of recalculation of photon energy with Doppler broadening
458
459 // *******************************************************
460 // | Determine ejected Compton electron direction |
461 // *******************************************************
462
463 // Calculate velocity of ejected Compton electron
464
465 G4double a_temp = eERecoil / electron_mass_c2;
466 G4double u_p_temp = sqrt(1 - (1 / (a_temp*a_temp)));
467
468 // Coefficients and terms from simulatenous equations
469 G4double sinAlpha = std::sin(e_alpha);
470 G4double cosAlpha = std::cos(e_alpha);
471 G4double sinBeta = std::sin(e_beta);
472 G4double cosBeta = std::cos(e_beta);
473
474 G4double gamma = 1.0 / sqrt(1 - (u_temp*u_temp));
475 G4double gamma_p = 1.0 / sqrt(1 - (u_p_temp*u_p_temp));
476
477 G4double var_A = pERecoil*u_p_temp*sinTheta;
478 G4double var_B = u_p_temp* (pERecoil*cosTheta-pEIncident);
479 G4double var_C = (pERecoil-pEIncident) - ( (pERecoil*pEIncident) / (gamma_p*electron_mass_c2))*(1 - cosTheta);
480
481 G4double var_D1 = gamma*electron_mass_c2*pERecoil;
482 G4double var_D2 = (1 - (u_temp*cosTheta*cosAlpha) - (u_temp*sinTheta*cosBeta*sinAlpha));
483 G4double var_D3 = ((electron_mass_c2*electron_mass_c2)*(gamma*gamma_p - 1)) - (gamma_p*electron_mass_c2*pERecoil);
484 G4double var_D = var_D1*var_D2 + var_D3;
485
486 G4double var_E1 = ((gamma*gamma_p)*(electron_mass_c2*electron_mass_c2)*(u_temp*u_p_temp)*cosAlpha);
487 G4double var_E2 = gamma_p*electron_mass_c2*pERecoil*u_p_temp*cosTheta;
488 G4double var_E = var_E1 - var_E2;
489
490 G4double var_F1 = ((gamma*gamma_p)*(electron_mass_c2*electron_mass_c2)*(u_temp*u_p_temp)*cosBeta*sinAlpha);
491 G4double var_F2 = (gamma_p*electron_mass_c2*pERecoil*u_p_temp*sinTheta);
492 G4double var_F = var_F1 - var_F2;
493
494 G4double var_G = (gamma*gamma_p)*(electron_mass_c2*electron_mass_c2)*(u_temp*u_p_temp)*sinBeta*sinAlpha;
495
496 // Two equations form a quadratic form of Wx^2 + Yx + Z = 0
497 // Coefficents and solution to quadratic
498 G4double var_W1 = (var_F*var_B - var_E*var_A)*(var_F*var_B - var_E*var_A);
499 G4double var_W2 = (var_G*var_G)*(var_A*var_A) + (var_G*var_G)*(var_B*var_B);
500 G4double var_W = var_W1 + var_W2;
501
502 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));
503
504 G4double var_Z1 = (var_A*var_D - var_F*var_C)*(var_A*var_D - var_F*var_C);
505 G4double var_Z2 = (var_G*var_G)*(var_C*var_C) - (var_G*var_G)*(var_A*var_A);
506 G4double var_Z = var_Z1 + var_Z2;
507 G4double diff1 = var_Y*var_Y;
508 G4double diff2 = 4*var_W*var_Z;
509 G4double diff = diff1 - diff2;
510
511 // Check if diff is less than zero, if so ensure it is due to FPE
512 //Determine number of digits (in decimal base) that G4double can accurately represent
513 G4double g4d_order = G4double(numeric_limits<G4double>::digits10);
514 G4double g4d_limit = std::pow(10.,-g4d_order);
515 //Confirm that diff less than zero is due FPE, i.e if abs of diff / diff1 and diff/ diff2 is less
516 //than 10^(-g4d_order), then set diff to zero
517
518 if ((diff < 0.0) && (abs(diff / diff1) < g4d_limit) && (abs(diff / diff2) < g4d_limit) )
519 {
520 diff = 0.0;
521 }
522
523
524 // Plus and minus of quadratic
525 G4double X_p = (-var_Y + sqrt (diff))/(2*var_W);
526 G4double X_m = (-var_Y - sqrt (diff))/(2*var_W);
527
528 // Floating point precision protection
529 // Check if X_p and X_m are greater than or less than 1 or -1, if so clean up FPE
530 // Issue due to propagation of FPE and only impacts 8th sig fig onwards
531 if(X_p >1){X_p=1;} if(X_p<-1){X_p=-1;}
532 if(X_m >1){X_m=1;} if(X_m<-1){X_m=-1;}
533 // End of FP protection
534
535 G4double ThetaE = 0.;
536 // Randomly sample one of the two possible solutions and determin theta angle of ejected Compton electron
537 G4double sol_select = G4UniformRand();
538
539 if (sol_select < 0.5)
540 {
541 ThetaE = std::acos(X_p);
542 }
543 if (sol_select > 0.5)
544 {
545 ThetaE = std::acos(X_m);
546 }
547 cosThetaE = std::cos(ThetaE);
548 sinThetaE = std::sin(ThetaE);
549 G4double Theta = std::acos(cosTheta);
550
551 //Calculate electron Phi
552 G4double iSinThetaE = std::sqrt(1+std::tan((pi/2.0)-ThetaE)*std::tan((pi/2.0)-ThetaE));
553 G4double iSinTheta = std::sqrt(1+std::tan((pi/2.0)-Theta)*std::tan((pi/2.0)-Theta));
554 G4double ivar_A = iSinTheta/ (pERecoil*u_p_temp);
555 // Trigs
556 cosPhiE = (var_C - var_B*cosThetaE)*(ivar_A*iSinThetaE);
557
558 // End of calculation of ejection Compton electron direction
559 //Fix for floating point errors
560
561 } while ( (iteration <= maxDopplerIterations) && (abs(cosPhiE) > 1));
562
563 // Revert to original if maximum number of iterations threshold has been reached
564 if (iteration >= maxDopplerIterations)
565 {
566 pERecoil = photonEnergy0 ;
567 bindingE = 0.;
568 dirx=0.0;
569 diry=0.0;
570 dirz=1.0;
571 }
572
573 // Set "scattered" photon direction and energy
574 G4ThreeVector photonDirection1(dirx,diry,dirz);
575 photonDirection1.rotateUz(photonDirection0);
576 fParticleChange->ProposeMomentumDirection(photonDirection1) ;
577
578 if (pERecoil > 0.)
579 {
580 fParticleChange->SetProposedKineticEnergy(pERecoil) ;
581
582 // Set ejected Compton electron direction and energy
583 G4double PhiE = std::acos(cosPhiE);
584 G4double eDirX = sinThetaE * std::cos(phi+PhiE);
585 G4double eDirY = sinThetaE * std::sin(phi+PhiE);
586 G4double eDirZ = cosThetaE;
587
588 G4double eKineticEnergy = pEIncident - pERecoil - bindingE;
589
590 G4ThreeVector eDirection(eDirX,eDirY,eDirZ);
591 eDirection.rotateUz(photonDirection0);
593 eDirection,eKineticEnergy) ;
594 fvect->push_back(dp);
595
596 }
597 else
598 {
599 fParticleChange->SetProposedKineticEnergy(0.);
600 fParticleChange->ProposeTrackStatus(fStopAndKill);
601 }
602
603 // sample deexcitation
604 //
605 if (verboseLevel > 3) {
606 G4cout << "Started atomic de-excitation " << fAtomDeexcitation << G4endl;
607 }
608
609 if(fAtomDeexcitation && iteration < maxDopplerIterations) {
610 G4int index = couple->GetIndex();
611 if(fAtomDeexcitation->CheckDeexcitationActiveRegion(index)) {
612 std::size_t nbefore = fvect->size();
614 const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
615 fAtomDeexcitation->GenerateParticles(fvect, shell, Z, index);
616 std::size_t nafter = fvect->size();
617 if(nafter > nbefore) {
618 for (std::size_t i=nbefore; i<nafter; ++i) {
619 //Check if there is enough residual energy
620 if (bindingE >= ((*fvect)[i])->GetKineticEnergy())
621 {
622 //Ok, this is a valid secondary: keep it
623 bindingE -= ((*fvect)[i])->GetKineticEnergy();
624 }
625 else
626 {
627 //Invalid secondary: not enough energy to create it!
628 //Keep its energy in the local deposit
629 delete (*fvect)[i];
630 (*fvect)[i]=nullptr;
631 }
632 }
633 }
634 }
635 }
636
637 //This should never happen
638 if(bindingE < 0.0)
639 G4Exception("G4LowEPComptonModel::SampleSecondaries()",
640 "em2051",FatalException,"Negative local energy deposit");
641
642 fParticleChange->ProposeLocalEnergyDeposit(bindingE);
643}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition G4Exp.hh:180
@ fStopAndKill
const G4double alpha2
#define G4UniformRand()
Definition Randomize.hh:52
G4double RandomSelectMomentum(G4int Z, G4int shellIndex) const
const G4ThreeVector & GetMomentumDirection() const
G4ParticleDefinition * GetDefinition() const
G4double GetKineticEnergy() const
static G4Electron * Electron()
Definition G4Electron.cc:91
G4double GetZ() const
Definition G4Element.hh:119
const G4String & GetName() const
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(const G4ThreeVector &Pfinal)
G4double BindingEnergy(G4int Z, G4int shellIndex) const
G4int SelectRandomShell(G4int Z) const
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)
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

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