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

#include <G4LowEPPolarizedComptonModel.hh>

+ Inheritance diagram for G4LowEPPolarizedComptonModel:

Public Member Functions

 G4LowEPPolarizedComptonModel (const G4ParticleDefinition *p=nullptr, const G4String &nam="LowEPComptonModel")
 
virtual ~G4LowEPPolarizedComptonModel ()
 
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
 
G4LowEPPolarizedComptonModeloperator= (const G4LowEPPolarizedComptonModel &right)=delete
 
 G4LowEPPolarizedComptonModel (const G4LowEPPolarizedComptonModel &)=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 91 of file G4LowEPPolarizedComptonModel.hh.

Constructor & Destructor Documentation

◆ G4LowEPPolarizedComptonModel() [1/2]

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

Definition at line 85 of file G4LowEPPolarizedComptonModel.cc.

87 : G4VEmModel(nam),isInitialised(false)
88{
89 verboseLevel=1 ;
90 // Verbosity scale:
91 // 0 = nothing
92 // 1 = warning for energy non-conservation
93 // 2 = details of energy budget
94 // 3 = calculation of cross sections, file openings, sampling of atoms
95 // 4 = entering in methods
96
97 if( verboseLevel>1 ) {
98 G4cout << "Low energy photon Compton model is constructed " << G4endl;
99 }
100
101 //Mark this model as "applicable" for atomic deexcitation
103
104 fParticleChange = nullptr;
105 fAtomDeexcitation = nullptr;
106}
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
void SetDeexcitationFlag(G4bool val)
G4VEmModel(const G4String &nam)
Definition G4VEmModel.cc:67

◆ ~G4LowEPPolarizedComptonModel()

G4LowEPPolarizedComptonModel::~G4LowEPPolarizedComptonModel ( )
virtual

Definition at line 110 of file G4LowEPPolarizedComptonModel.cc.

111{
112 if(IsMaster()) {
113 delete shellData;
114 shellData = nullptr;
115 delete profileData;
116 profileData = nullptr;
117 }
118}
G4bool IsMaster() const

◆ G4LowEPPolarizedComptonModel() [2/2]

G4LowEPPolarizedComptonModel::G4LowEPPolarizedComptonModel ( const G4LowEPPolarizedComptonModel & )
delete

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 244 of file G4LowEPPolarizedComptonModel.cc.

248{
249 if (verboseLevel > 3) {
250 G4cout << "G4LowEPPolarizedComptonModel::ComputeCrossSectionPerAtom()"
251 << G4endl;
252 }
253 G4double cs = 0.0;
254
255 if (GammaEnergy < LowEnergyLimit()) { return 0.0; }
256
257 G4int intZ = G4lrint(Z);
258 if(intZ < 1 || intZ > maxZ) { return cs; }
259
260 G4PhysicsFreeVector* pv = data[intZ];
261
262 // if element was not initialised
263 // do initialisation safely for MT mode
264 if(!pv)
265 {
266 InitialiseForElement(0, intZ);
267 pv = data[intZ];
268 if(!pv) { return cs; }
269 }
270
271 G4int n = G4int(pv->GetVectorLength() - 1);
272 G4double e1 = pv->Energy(0);
273 G4double e2 = pv->Energy(n);
274
275 if(GammaEnergy <= e1) { cs = GammaEnergy/(e1*e1)*pv->Value(e1); }
276 else if(GammaEnergy <= e2) { cs = pv->Value(GammaEnergy)/GammaEnergy; }
277 else if(GammaEnergy > e2) { cs = pv->Value(e2)/GammaEnergy; }
278
279 return cs;
280}
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 G4LowEPPolarizedComptonModel::Initialise ( const G4ParticleDefinition * particle,
const G4DataVector & cuts )
overridevirtual

Implements G4VEmModel.

Definition at line 122 of file G4LowEPPolarizedComptonModel.cc.

124{
125 if (verboseLevel > 1) {
126 G4cout << "Calling G4LowEPPolarizedComptonModel::Initialise()" << G4endl;
127 }
128
129 // Initialise element selector
130 if(IsMaster()) {
131 // Access to elements
132 const char* path = G4FindDataDir("G4LEDATA");
133
134 G4ProductionCutsTable* theCoupleTable =
136 G4int numOfCouples = (G4int)theCoupleTable->GetTableSize();
137
138 for(G4int i=0; i<numOfCouples; ++i) {
139 const G4Material* material =
140 theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
141 const G4ElementVector* theElementVector = material->GetElementVector();
142 std::size_t nelm = material->GetNumberOfElements();
143
144 for (std::size_t j=0; j<nelm; ++j) {
145 G4int Z = G4lrint((*theElementVector)[j]->GetZ());
146 if(Z < 1) { Z = 1; }
147 else if(Z > maxZ){ Z = maxZ; }
148
149 if( (!data[Z]) ) { ReadData(Z, path); }
150 }
151 }
152
153 // For Doppler broadening
154 if(!shellData) {
155 shellData = new G4ShellData();
156 shellData->SetOccupancyData();
157 G4String file = "/doppler/shell-doppler";
158 shellData->LoadData(file);
159 }
160 if(!profileData) { profileData = new G4DopplerProfile(); }
161
162 InitialiseElementSelectors(particle, cuts);
163 }
164
165 if (verboseLevel > 2) {
166 G4cout << "Loaded cross section files" << G4endl;
167 }
168
169 if( verboseLevel>1 ) {
170 G4cout << "G4LowEPPolarizedComptonModel is initialized " << G4endl
171 << "Energy range: "
172 << LowEnergyLimit() / eV << " eV - "
173 << HighEnergyLimit() / GeV << " GeV"
174 << G4endl;
175 }
176
177 if(isInitialised) { return; }
178
179 fParticleChange = GetParticleChangeForGamma();
180 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
181 isInitialised = true;
182}
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 G4LowEPPolarizedComptonModel::InitialiseForElement ( const G4ParticleDefinition * ,
G4int Z )
overridevirtual

Reimplemented from G4VEmModel.

Definition at line 684 of file G4LowEPPolarizedComptonModel.cc.

686{
687 G4AutoLock l(&LowEPPolarizedComptonModelMutex);
688 if(!data[Z]) { ReadData(Z); }
689 l.unlock();
690}

Referenced by ComputeCrossSectionPerAtom().

◆ InitialiseLocal()

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

Reimplemented from G4VEmModel.

Definition at line 186 of file G4LowEPPolarizedComptonModel.cc.

188{
190}
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
std::vector< G4EmElementSelector * > * GetElementSelectors()

◆ operator=()

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

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 284 of file G4LowEPPolarizedComptonModel.cc.

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