Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions
G4LivermoreComptonModel Class Reference

#include <G4LivermoreComptonModel.hh>

Inheritance diagram for G4LivermoreComptonModel:
G4VEmModel

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 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 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 *)
 
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)
 
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)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
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
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
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 SetForceBuildTable (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

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 std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Detailed Description

Definition at line 48 of file G4LivermoreComptonModel.hh.

Constructor & Destructor Documentation

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

Definition at line 69 of file G4LivermoreComptonModel.cc.

References python.hepunit::eV, G4cout, G4endl, and G4VEmModel::SetDeexcitationFlag().

71  : G4VEmModel(nam),isInitialised(false)
72 {
73  lowestEnergy = 10 * eV;
74 
75  verboseLevel=1 ;
76  // Verbosity scale:
77  // 0 = nothing
78  // 1 = warning for energy non-conservation
79  // 2 = details of energy budget
80  // 3 = calculation of cross sections, file openings, sampling of atoms
81  // 4 = entering in methods
82 
83  if( verboseLevel>1 ) {
84  G4cout << "Livermore Compton model is constructed " << G4endl;
85  }
86 
87  //Mark this model as "applicable" for atomic deexcitation
88  SetDeexcitationFlag(true);
89 
90  fParticleChange = 0;
91  fAtomDeexcitation = 0;
92 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:65
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
void SetDeexcitationFlag(G4bool val)
Definition: G4VEmModel.hh:739
G4LivermoreComptonModel::~G4LivermoreComptonModel ( )
virtual

Definition at line 96 of file G4LivermoreComptonModel.cc.

References G4VEmModel::IsMaster().

97 {
98  if(IsMaster()) {
99  delete shellData;
100  shellData = 0;
101  delete profileData;
102  profileData = 0;
103  }
104 }
G4bool IsMaster() const
Definition: G4VEmModel.hh:676

Member Function Documentation

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 236 of file G4LivermoreComptonModel.cc.

References G4PhysicsVector::Energy(), G4cout, G4endl, G4lrint(), G4PhysicsVector::GetVectorLength(), InitialiseForElement(), n, and G4PhysicsVector::Value().

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

Implements G4VEmModel.

Definition at line 108 of file G4LivermoreComptonModel.cc.

References G4LossTableManager::AtomDeexcitation(), python.hepunit::eV, G4cout, G4endl, G4lrint(), G4Material::GetElementVector(), G4MaterialCutsCouple::GetMaterial(), G4ProductionCutsTable::GetMaterialCutsCouple(), G4Material::GetNumberOfElements(), G4VEmModel::GetParticleChangeForGamma(), G4ProductionCutsTable::GetProductionCutsTable(), G4ProductionCutsTable::GetTableSize(), python.hepunit::GeV, G4VEmModel::HighEnergyLimit(), G4VEmModel::InitialiseElementSelectors(), G4LossTableManager::Instance(), G4VEmModel::IsMaster(), G4ShellData::LoadData(), G4VEmModel::LowEnergyLimit(), eplot::material, and G4ShellData::SetOccupancyData().

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

Reimplemented from G4VEmModel.

Definition at line 549 of file G4LivermoreComptonModel.cc.

References G4TemplateAutoLock< M, L, U >::unlock().

Referenced by ComputeCrossSectionPerAtom().

551 {
552  G4AutoLock l(&LivermoreComptonModelMutex);
553  // G4cout << "G4LivermoreComptonModel::InitialiseForElement Z= "
554  // << Z << G4endl;
555  if(!data[Z]) { ReadData(Z); }
556  l.unlock();
557 }
const XML_Char const XML_Char * data
void G4LivermoreComptonModel::InitialiseLocal ( const G4ParticleDefinition ,
G4VEmModel masterModel 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 175 of file G4LivermoreComptonModel.cc.

References G4VEmModel::GetElementSelectors(), and G4VEmModel::SetElementSelectors().

177 {
179 }
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:760
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:768
void G4LivermoreComptonModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle aDynamicGamma,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 277 of file G4LivermoreComptonModel.cc.

References G4ShellData::BindingEnergy(), python.hepunit::c_light, G4VAtomDeexcitation::CheckDeexcitationActiveRegion(), python.hepunit::cm, G4Electron::Electron(), python.hepunit::electron_mass_c2, python.hepunit::fine_structure_const, fStopAndKill, G4cout, G4endl, G4Exp(), G4Log(), G4lrint(), G4UniformRand, G4VAtomDeexcitation::GenerateParticles(), G4VAtomDeexcitation::GetAtomicShell(), G4DynamicParticle::GetDefinition(), G4MaterialCutsCouple::GetIndex(), G4DynamicParticle::GetKineticEnergy(), G4MaterialCutsCouple::GetMaterial(), G4DynamicParticle::GetMomentumDirection(), G4Material::GetName(), G4Element::GetZ(), python.hepunit::h_Planck, python.hepunit::keV, python.hepunit::MeV, G4VParticleChange::ProposeLocalEnergyDeposit(), G4ParticleChangeForGamma::ProposeMomentumDirection(), G4VParticleChange::ProposeTrackStatus(), G4DopplerProfile::RandomSelectMomentum(), CLHEP::Hep3Vector::rotateUz(), G4VEmModel::SelectRandomAtom(), G4ShellData::SelectRandomShell(), G4ParticleChangeForGamma::SetProposedKineticEnergy(), python.hepunit::twopi, var(), and test::x.

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

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