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

#include <G4LivermoreComptonModel.hh>

+ Inheritance diagram for G4LivermoreComptonModel:

Public Member Functions

 G4LivermoreComptonModel (const G4ParticleDefinition *p=0, const G4String &nam="LivermoreCompton")
 
virtual ~G4LivermoreComptonModel ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual void InitialiseForElement (const G4ParticleDefinition *, G4int Z)
 
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 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 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 46 of file G4LivermoreComptonModel.hh.

Constructor & Destructor Documentation

◆ G4LivermoreComptonModel()

G4LivermoreComptonModel::G4LivermoreComptonModel ( const G4ParticleDefinition p = 0,
const G4String nam = "LivermoreCompton" 
)

Definition at line 67 of file G4LivermoreComptonModel.cc.

69 : G4VEmModel(nam),isInitialised(false)
70{
71 verboseLevel=1 ;
72 // Verbosity scale:
73 // 0 = nothing
74 // 1 = warning for energy non-conservation
75 // 2 = details of energy budget
76 // 3 = calculation of cross sections, file openings, sampling of atoms
77 // 4 = entering in methods
78
79 if( verboseLevel>1 ) {
80 G4cout << "Livermore Compton model is constructed " << G4endl;
81 }
82
83 //Mark this model as "applicable" for atomic deexcitation
85
86 fParticleChange = 0;
87 fAtomDeexcitation = 0;
88}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:813

◆ ~G4LivermoreComptonModel()

G4LivermoreComptonModel::~G4LivermoreComptonModel ( )
virtual

Definition at line 92 of file G4LivermoreComptonModel.cc.

93{
94 if(IsMaster()) {
95 delete shellData;
96 shellData = 0;
97 delete profileData;
98 profileData = 0;
99 for(G4int i=0; i<maxZ; ++i) {
100 if(data[i]) {
101 delete data[i];
102 data[i] = 0;
103 }
104 }
105 }
106}
int G4int
Definition: G4Types.hh:85
G4bool IsMaster() const
Definition: G4VEmModel.hh:736

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 238 of file G4LivermoreComptonModel.cc.

242{
243 if (verboseLevel > 3) {
244 G4cout << "G4LivermoreComptonModel::ComputeCrossSectionPerAtom()"
245 << G4endl;
246 }
247 G4double cs = 0.0;
248
249 if (GammaEnergy < LowEnergyLimit()) { return 0.0; }
250
251 G4int intZ = G4lrint(Z);
252 if(intZ < 1 || intZ > maxZ) { return cs; }
253
254 G4LPhysicsFreeVector* pv = data[intZ];
255
256 // if element was not initialised
257 // do initialisation safely for MT mode
258 if(!pv)
259 {
260 InitialiseForElement(0, intZ);
261 pv = data[intZ];
262 if(!pv) { return cs; }
263 }
264
265 G4int n = pv->GetVectorLength() - 1;
266 G4double e1 = pv->Energy(0);
267 G4double e2 = pv->Energy(n);
268
269 if(GammaEnergy <= e1) { cs = GammaEnergy/(e1*e1)*pv->Value(e1); }
270 else if(GammaEnergy <= e2) { cs = pv->Value(GammaEnergy)/GammaEnergy; }
271 else if(GammaEnergy > e2) { cs = pv->Value(e2)/GammaEnergy; }
272
273 return cs;
274}
double G4double
Definition: G4Types.hh:83
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
G4double Energy(std::size_t index) const
G4double Value(G4double theEnergy, std::size_t &lastidx) const
std::size_t GetVectorLength() const
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
int G4lrint(double ad)
Definition: templates.hh:134

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 110 of file G4LivermoreComptonModel.cc.

112{
113 if (verboseLevel > 1) {
114 G4cout << "Calling G4LivermoreComptonModel::Initialise()" << G4endl;
115 }
116
117 // Initialise element selector
118
119 if(IsMaster()) {
120
121 // Access to elements
122
123 char* path = std::getenv("G4LEDATA");
124
125 G4ProductionCutsTable* theCoupleTable =
127 G4int numOfCouples = theCoupleTable->GetTableSize();
128
129 for(G4int i=0; i<numOfCouples; ++i) {
130 const G4Material* material =
131 theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
132 const G4ElementVector* theElementVector = material->GetElementVector();
133 G4int nelm = material->GetNumberOfElements();
134
135 for (G4int j=0; j<nelm; ++j) {
136 G4int Z = G4lrint((*theElementVector)[j]->GetZ());
137 if(Z < 1) { Z = 1; }
138 else if(Z > maxZ){ Z = maxZ; }
139
140 if( (!data[Z]) ) { ReadData(Z, path); }
141 }
142 }
143
144 // For Doppler broadening
145 if(!shellData) {
146 shellData = new G4ShellData();
147 shellData->SetOccupancyData();
148 G4String file = "/doppler/shell-doppler";
149 shellData->LoadData(file);
150 }
151 if(!profileData) { profileData = new G4DopplerProfile(); }
152
153 InitialiseElementSelectors(particle, cuts);
154 }
155
156 if (verboseLevel > 2) {
157 G4cout << "Loaded cross section files" << G4endl;
158 }
159
160 if( verboseLevel>1 ) {
161 G4cout << "G4LivermoreComptonModel is initialized " << G4endl
162 << "Energy range: "
163 << LowEnergyLimit() / eV << " eV - "
164 << HighEnergyLimit() / GeV << " GeV"
165 << G4endl;
166 }
167 //
168 if(isInitialised) { return; }
169
170 fParticleChange = GetParticleChangeForGamma();
171 fAtomDeexcitation = G4LossTableManager::Instance()->AtomDeexcitation();
172 isInitialised = true;
173}
std::vector< G4Element * > G4ElementVector
static G4LossTableManager * Instance()
G4VAtomDeexcitation * AtomDeexcitation()
const G4Material * GetMaterial() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:188
size_t GetNumberOfElements() const
Definition: G4Material.hh:184
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
void SetOccupancyData()
Definition: G4ShellData.hh:69
void LoadData(const G4String &fileName)
Definition: G4ShellData.cc:233
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148

◆ InitialiseForElement()

void G4LivermoreComptonModel::InitialiseForElement ( const G4ParticleDefinition ,
G4int  Z 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 563 of file G4LivermoreComptonModel.cc.

565{
566 G4AutoLock l(&LivermoreComptonModelMutex);
567 // G4cout << "G4LivermoreComptonModel::InitialiseForElement Z= "
568 // << Z << G4endl;
569 if(!data[Z]) { ReadData(Z); }
570 l.unlock();
571}

Referenced by ComputeCrossSectionPerAtom().

◆ InitialiseLocal()

void G4LivermoreComptonModel::InitialiseLocal ( const G4ParticleDefinition ,
G4VEmModel masterModel 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 177 of file G4LivermoreComptonModel.cc.

179{
181}
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:842
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:834

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 279 of file G4LivermoreComptonModel.cc.

284{
285
286 // The scattered gamma energy is sampled according to Klein - Nishina
287 // formula then accepted or rejected depending on the Scattering Function
288 // multiplied by factor from Klein - Nishina formula.
289 // Expression of the angular distribution as Klein Nishina
290 // angular and energy distribution and Scattering fuctions is taken from
291 // D. E. Cullen "A simple model of photon transport" Nucl. Instr. Meth.
292 // Phys. Res. B 101 (1995). Method of sampling with form factors is different
293 // data are interpolated while in the article they are fitted.
294 // Reference to the article is from J. Stepanek New Photon, Positron
295 // and Electron Interaction Data for GEANT in Energy Range from 1 eV to 10
296 // TeV (draft).
297 // The random number techniques of Butcher & Messel are used
298 // (Nucl Phys 20(1960),15).
299
300 G4double photonEnergy0 = aDynamicGamma->GetKineticEnergy();
301
302 if (verboseLevel > 3) {
303 G4cout << "G4LivermoreComptonModel::SampleSecondaries() E(MeV)= "
304 << photonEnergy0/MeV << " in " << couple->GetMaterial()->GetName()
305 << G4endl;
306 }
307
308 // do nothing below the threshold
309 // should never get here because the XS is zero below the limit
310 if (photonEnergy0 < LowEnergyLimit())
311 return ;
312
313 G4double e0m = photonEnergy0 / electron_mass_c2 ;
314 G4ParticleMomentum photonDirection0 = aDynamicGamma->GetMomentumDirection();
315
316 // Select randomly one element in the current material
317 const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
318 const G4Element* elm = SelectRandomAtom(couple,particle,photonEnergy0);
319
320 G4int Z = G4lrint(elm->GetZ());
321
322 G4double epsilon0Local = 1. / (1. + 2. * e0m);
323 G4double epsilon0Sq = epsilon0Local * epsilon0Local;
324 G4double alpha1 = -G4Log(epsilon0Local);
325 G4double alpha2 = 0.5 * (1. - epsilon0Sq);
326
327 G4double wlPhoton = h_Planck*c_light/photonEnergy0;
328
329 // Sample the energy of the scattered photon
331 G4double epsilonSq;
332 G4double oneCosT;
333 G4double sinT2;
334 G4double gReject;
335
336 if (verboseLevel > 3) {
337 G4cout << "Started loop to sample gamma energy" << G4endl;
338 }
339
340 do {
341 if ( alpha1/(alpha1+alpha2) > G4UniformRand())
342 {
343 epsilon = G4Exp(-alpha1 * G4UniformRand());
344 epsilonSq = epsilon * epsilon;
345 }
346 else
347 {
348 epsilonSq = epsilon0Sq + (1. - epsilon0Sq) * G4UniformRand();
349 epsilon = std::sqrt(epsilonSq);
350 }
351
352 oneCosT = (1. - epsilon) / ( epsilon * e0m);
353 sinT2 = oneCosT * (2. - oneCosT);
354 G4double x = std::sqrt(oneCosT/2.) * cm / wlPhoton;
355 G4double scatteringFunction = ComputeScatteringFunction(x, Z);
356 gReject = (1. - epsilon * sinT2 / (1. + epsilonSq)) * scatteringFunction;
357
358 } while(gReject < G4UniformRand()*Z);
359
360 G4double cosTheta = 1. - oneCosT;
361 G4double sinTheta = std::sqrt (sinT2);
362 G4double phi = twopi * G4UniformRand() ;
363 G4double dirx = sinTheta * std::cos(phi);
364 G4double diry = sinTheta * std::sin(phi);
365 G4double dirz = cosTheta ;
366
367 // Doppler broadening - Method based on:
368 // Y. Namito, S. Ban and H. Hirayama,
369 // "Implementation of the Doppler Broadening of a Compton-Scattered Photon
370 // into the EGS4 Code", NIM A 349, pp. 489-494, 1994
371
372 // Maximum number of sampling iterations
373 static G4int maxDopplerIterations = 1000;
374 G4double bindingE = 0.;
375 G4double photonEoriginal = epsilon * photonEnergy0;
376 G4double photonE = -1.;
377 G4int iteration = 0;
378 G4double eMax = photonEnergy0;
379
380 G4int shellIdx = 0;
381
382 if (verboseLevel > 3) {
383 G4cout << "Started loop to sample broading" << G4endl;
384 }
385
386 do
387 {
388 ++iteration;
389 // Select shell based on shell occupancy
390 shellIdx = shellData->SelectRandomShell(Z);
391 bindingE = shellData->BindingEnergy(Z,shellIdx);
392
393 if (verboseLevel > 3) {
394 G4cout << "Shell ID= " << shellIdx
395 << " Ebind(keV)= " << bindingE/keV << G4endl;
396 }
397
398 eMax = photonEnergy0 - bindingE;
399
400 // Randomly sample bound electron momentum
401 // (memento: the data set is in Atomic Units)
402 G4double pSample = profileData->RandomSelectMomentum(Z,shellIdx);
403 if (verboseLevel > 3) {
404 G4cout << "pSample= " << pSample << G4endl;
405 }
406 // Rescale from atomic units
407 G4double pDoppler = pSample * fine_structure_const;
408 G4double pDoppler2 = pDoppler * pDoppler;
409 G4double var2 = 1. + oneCosT * e0m;
410 G4double var3 = var2*var2 - pDoppler2;
411 G4double var4 = var2 - pDoppler2 * cosTheta;
412 G4double var = var4*var4 - var3 + pDoppler2 * var3;
413 if (var > 0.)
414 {
415 G4double varSqrt = std::sqrt(var);
416 G4double scale = photonEnergy0 / var3;
417 // Random select either root
418 if (G4UniformRand() < 0.5) { photonE = (var4 - varSqrt) * scale; }
419 else { photonE = (var4 + varSqrt) * scale; }
420 }
421 else
422 {
423 photonE = -1.;
424 }
425 } while (iteration <= maxDopplerIterations && photonE > eMax);
426 // (photonE < 0. || photonE > eMax || photonE < eMax*G4UniformRand()) );
427
428
429 // End of recalculation of photon energy with Doppler broadening
430 // Revert to original if maximum number of iterations threshold
431 // has been reached
432
433 if (iteration >= maxDopplerIterations)
434 {
435 photonE = photonEoriginal;
436 bindingE = 0.;
437 }
438
439 // Update G4VParticleChange for the scattered photon
440
441 G4ThreeVector photonDirection1(dirx,diry,dirz);
442 photonDirection1.rotateUz(photonDirection0);
443 fParticleChange->ProposeMomentumDirection(photonDirection1) ;
444
445 G4double photonEnergy1 = photonE;
446
447 if (photonEnergy1 > 0.) {
448 fParticleChange->SetProposedKineticEnergy(photonEnergy1) ;
449
450 } else {
451 // photon absorbed
452 photonEnergy1 = 0.;
453 fParticleChange->SetProposedKineticEnergy(0.) ;
454 fParticleChange->ProposeTrackStatus(fStopAndKill);
455 fParticleChange->ProposeLocalEnergyDeposit(photonEnergy0);
456 return;
457 }
458
459 // Kinematics of the scattered electron
460 G4double eKineticEnergy = photonEnergy0 - photonEnergy1 - bindingE;
461
462 // protection against negative final energy: no e- is created
463 if(eKineticEnergy < 0.0) {
464 fParticleChange->ProposeLocalEnergyDeposit(photonEnergy0 - photonEnergy1);
465 return;
466 }
467
468 G4double eTotalEnergy = eKineticEnergy + electron_mass_c2;
469
470 G4double electronE = photonEnergy0 * (1. - epsilon) + electron_mass_c2;
471 G4double electronP2 =
472 electronE*electronE - electron_mass_c2*electron_mass_c2;
473 G4double sinThetaE = -1.;
474 G4double cosThetaE = 0.;
475 if (electronP2 > 0.)
476 {
477 cosThetaE = (eTotalEnergy + photonEnergy1 )*
478 (1. - epsilon) / std::sqrt(electronP2);
479 sinThetaE = -1. * sqrt(1. - cosThetaE * cosThetaE);
480 }
481
482 G4double eDirX = sinThetaE * std::cos(phi);
483 G4double eDirY = sinThetaE * std::sin(phi);
484 G4double eDirZ = cosThetaE;
485
486 G4ThreeVector eDirection(eDirX,eDirY,eDirZ);
487 eDirection.rotateUz(photonDirection0);
489 eDirection,eKineticEnergy) ;
490 fvect->push_back(dp);
491
492 // sample deexcitation
493 //
494
495 if (verboseLevel > 3) {
496 G4cout << "Started atomic de-excitation " << fAtomDeexcitation << G4endl;
497 }
498
499 if(fAtomDeexcitation && iteration < maxDopplerIterations) {
500 G4int index = couple->GetIndex();
501 if(fAtomDeexcitation->CheckDeexcitationActiveRegion(index)) {
502 size_t nbefore = fvect->size();
504 const G4AtomicShell* shell = fAtomDeexcitation->GetAtomicShell(Z, as);
505 fAtomDeexcitation->GenerateParticles(fvect, shell, Z, index);
506 size_t nafter = fvect->size();
507 if(nafter > nbefore) {
508 for (size_t i=nbefore; i<nafter; ++i) {
509 //Check if there is enough residual energy
510 if (bindingE >= ((*fvect)[i])->GetKineticEnergy())
511 {
512 //Ok, this is a valid secondary: keep it
513 bindingE -= ((*fvect)[i])->GetKineticEnergy();
514 }
515 else
516 {
517 //Invalid secondary: not enough energy to create it!
518 //Keep its energy in the local deposit
519 delete (*fvect)[i];
520 (*fvect)[i]=0;
521 }
522 }
523 }
524 }
525 }
526 //This should never happen
527 if(bindingE < 0.0)
528 G4Exception("G4LivermoreComptonModel::SampleSecondaries()",
529 "em2050",FatalException,"Negative local energy deposit");
530
531 fParticleChange->ProposeLocalEnergyDeposit(bindingE);
532}
G4AtomicShellEnumerator
double epsilon(double density, double temperature)
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:179
G4double G4Log(G4double x)
Definition: G4Log.hh:226
@ fStopAndKill
#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:93
G4double GetZ() const
Definition: G4Element.hh:130
const G4String & GetName() const
Definition: G4Material.hh:175
void SetProposedKineticEnergy(G4double proposedKinEnergy)
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
G4double BindingEnergy(G4int Z, G4int shellIndex) const
Definition: G4ShellData.cc:165
G4int SelectRandomShell(G4int Z) const
Definition: G4ShellData.cc:362
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)
Definition: G4VEmModel.hh:570
void ProposeTrackStatus(G4TrackStatus status)
void ProposeLocalEnergyDeposit(G4double anEnergyPart)

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