Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4LivermorePolarizedGammaConversionModel Class Reference

#include <G4LivermorePolarizedGammaConversionModel.hh>

Inheritance diagram for G4LivermorePolarizedGammaConversionModel:
Collaboration diagram for G4LivermorePolarizedGammaConversionModel:

Public Member Functions

 G4LivermorePolarizedGammaConversionModel (const G4ParticleDefinition *p=0, const G4String &nam="LivermorePolarizedGammaConversion")
 
virtual ~G4LivermorePolarizedGammaConversionModel ()
 
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)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double)
 
- 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 *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
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)
 
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=nullptr)
 
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
 
const G4IsotopeGetCurrentIsotope () const
 
G4bool IsLocked () const
 
void SetLocked (G4bool)
 

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
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 43 of file G4LivermorePolarizedGammaConversionModel.hh.

Constructor & Destructor Documentation

G4LivermorePolarizedGammaConversionModel::G4LivermorePolarizedGammaConversionModel ( const G4ParticleDefinition p = 0,
const G4String nam = "LivermorePolarizedGammaConversion" 
)

Definition at line 51 of file G4LivermorePolarizedGammaConversionModel.cc.

53  :G4VEmModel(nam), isInitialised(false),smallEnergy(2.*MeV)
54 {
55  fParticleChange = nullptr;
56  lowEnergyLimit = 2*electron_mass_c2;
57 
58  Phi=0.;
59  Psi=0.;
60 
61  verboseLevel= 0;
62  // Verbosity scale:
63  // 0 = nothing
64  // 1 = calculation of cross sections, file openings, samping of atoms
65  // 2 = entering in methods
66 
67  if(verboseLevel > 0) {
68  G4cout << "Livermore Polarized GammaConversion is constructed "
69  << G4endl;
70  }
71 
72 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
static constexpr double electron_mass_c2
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
G4LivermorePolarizedGammaConversionModel::~G4LivermorePolarizedGammaConversionModel ( )
virtual

Definition at line 76 of file G4LivermorePolarizedGammaConversionModel.cc.

77 {
78  if(IsMaster()) {
79  for(G4int i=0; i<maxZ; ++i) {
80  if(data[i]) {
81  delete data[i];
82  data[i] = 0;
83  }
84  }
85  }
86 }
G4bool IsMaster() const
Definition: G4VEmModel.hh:717
int G4int
Definition: G4Types.hh:78
const XML_Char const XML_Char * data
Definition: expat.h:268

Here is the call graph for this function:

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 221 of file G4LivermorePolarizedGammaConversionModel.cc.

226 {
227  if (verboseLevel > 1) {
228  G4cout << "G4LivermorePolarizedGammaConversionModel::ComputeCrossSectionPerAtom()"
229  << G4endl;
230  }
231  if (GammaEnergy < lowEnergyLimit) { return 0.0; }
232 
233  G4double xs = 0.0;
234 
235  G4int intZ=G4int(Z);
236 
237  if(intZ < 1 || intZ > maxZ) { return xs; }
238 
239  G4LPhysicsFreeVector* pv = data[intZ];
240 
241  // if element was not initialised
242  // do initialisation safely for MT mode
243  if(!pv)
244  {
245  InitialiseForElement(0, intZ);
246  pv = data[intZ];
247  if(!pv) { return xs; }
248  }
249  // x-section is taken from the table
250  xs = pv->Value(GammaEnergy);
251 
252  if(verboseLevel > 0)
253  {
254  G4int n = pv->GetVectorLength() - 1;
255  G4cout << "****** DEBUG: tcs value for Z=" << Z << " at energy (MeV)="
256  << GammaEnergy/MeV << G4endl;
257  G4cout << " cs (Geant4 internal unit)=" << xs << G4endl;
258  G4cout << " -> first cs value in EADL data file (iu) =" << (*pv)[0] << G4endl;
259  G4cout << " -> last cs value in EADL data file (iu) =" << (*pv)[n] << G4endl;
260  G4cout << "*********************************************************" << G4endl;
261  }
262 
263  return xs;
264 }
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)
size_t GetVectorLength() const
int G4int
Definition: G4Types.hh:78
const XML_Char const XML_Char * data
Definition: expat.h:268
G4GLOB_DLL std::ostream G4cout
G4double Value(G4double theEnergy, size_t &lastidx) const
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 90 of file G4LivermorePolarizedGammaConversionModel.cc.

92 {
93  if (verboseLevel > 1)
94  {
95  G4cout << "Calling1 G4LivermorePolarizedGammaConversionModel::Initialise()"
96  << G4endl
97  << "Energy range: "
98  << LowEnergyLimit() / MeV << " MeV - "
99  << HighEnergyLimit() / GeV << " GeV"
100  << G4endl;
101  }
102 
103 
104  if(IsMaster())
105  {
106 
107  // Initialise element selector
108 
109  InitialiseElementSelectors(particle, cuts);
110 
111  // Access to elements
112 
113  char* path = getenv("G4LEDATA");
114 
115  G4ProductionCutsTable* theCoupleTable =
117 
118  G4int numOfCouples = theCoupleTable->GetTableSize();
119 
120  for(G4int i=0; i<numOfCouples; ++i)
121  {
122  const G4Material* material =
123  theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
124  const G4ElementVector* theElementVector = material->GetElementVector();
125  G4int nelm = material->GetNumberOfElements();
126 
127  for (G4int j=0; j<nelm; ++j)
128  {
129  G4int Z = (G4int)(*theElementVector)[j]->GetZ();
130  if(Z < 1) { Z = 1; }
131  else if(Z > maxZ) { Z = maxZ; }
132  if(!data[Z]) { ReadData(Z, path); }
133  }
134  }
135  }
136  if(isInitialised) { return; }
137  fParticleChange = GetParticleChangeForGamma();
138  isInitialised = true;
139 
140 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:640
std::vector< G4Element * > G4ElementVector
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:146
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:633
G4bool IsMaster() const
Definition: G4VEmModel.hh:717
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:190
int G4int
Definition: G4Types.hh:78
const XML_Char const XML_Char * data
Definition: expat.h:268
G4GLOB_DLL std::ostream G4cout
static G4ProductionCutsTable * GetProductionCutsTable()
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
static constexpr double GeV
Definition: G4SIunits.hh:217
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
size_t GetNumberOfElements() const
Definition: G4Material.hh:186
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:132
const G4Material * GetMaterial() const

Here is the call graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 1087 of file G4LivermorePolarizedGammaConversionModel.cc.

1090 {
1091  G4AutoLock l(&LivermorePolarizedGammaConversionModelMutex);
1092  // G4cout << "G4LivermorePolarizedGammaConversionModel::InitialiseForElement Z= "
1093  // << Z << G4endl;
1094  if(!data[Z]) { ReadData(Z); }
1095  l.unlock();
1096 }
const XML_Char const XML_Char * data
Definition: expat.h:268

Here is the call graph for this function:

Here is the caller graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 145 of file G4LivermorePolarizedGammaConversionModel.cc.

147 {
149 }
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:801
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:809

Here is the call graph for this function:

G4double G4LivermorePolarizedGammaConversionModel::MinPrimaryEnergy ( const G4Material ,
const G4ParticleDefinition ,
G4double   
)
virtual

Reimplemented from G4VEmModel.

Definition at line 153 of file G4LivermorePolarizedGammaConversionModel.cc.

155 {
156  return lowEnergyLimit;
157 }
void G4LivermorePolarizedGammaConversionModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle aDynamicGamma,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 269 of file G4LivermorePolarizedGammaConversionModel.cc.

274 {
275 
276  // Fluorescence generated according to:
277  // J. Stepanek ,"A program to determine the radiation spectra due to a single atomic
278  // subshell ionisation by a particle or due to deexcitation or decay of radionuclides",
279  // Comp. Phys. Comm. 1206 pp 1-1-9 (1997)
280 
281  if (verboseLevel > 3)
282  G4cout << "Calling SampleSecondaries() of G4LivermorePolarizedGammaConversionModel" << G4endl;
283 
284  G4double photonEnergy = aDynamicGamma->GetKineticEnergy();
285  // Within energy limit?
286 
287  if(photonEnergy <= lowEnergyLimit)
288  {
289  fParticleChange->ProposeTrackStatus(fStopAndKill);
290  fParticleChange->SetProposedKineticEnergy(0.);
291  return;
292  }
293 
294 
295  G4ThreeVector gammaPolarization0 = aDynamicGamma->GetPolarization();
296  G4ThreeVector gammaDirection0 = aDynamicGamma->GetMomentumDirection();
297 
298  // Make sure that the polarization vector is perpendicular to the
299  // gamma direction. If not
300 
301  if(!(gammaPolarization0.isOrthogonal(gammaDirection0, 1e-6))||(gammaPolarization0.mag()==0))
302  { // only for testing now
303  gammaPolarization0 = GetRandomPolarization(gammaDirection0);
304  }
305  else
306  {
307  if ( gammaPolarization0.howOrthogonal(gammaDirection0) != 0)
308  {
309  gammaPolarization0 = GetPerpendicularPolarization(gammaDirection0, gammaPolarization0);
310  }
311  }
312 
313  // End of Protection
314 
315 
316  G4double epsilon ;
317  G4double epsilon0Local = electron_mass_c2 / photonEnergy ;
318 
319  // Do it fast if photon energy < 2. MeV
320 
321  if (photonEnergy < smallEnergy )
322  {
323  epsilon = epsilon0Local + (0.5 - epsilon0Local) * G4UniformRand();
324  }
325  else
326  {
327 
328  // Select randomly one element in the current material
329 
330  const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
331  const G4Element* element = SelectRandomAtom(couple,particle,photonEnergy);
332 
333 
334  if (element == 0)
335  {
336  G4cout << "G4LivermorePolarizedGammaConversionModel::SampleSecondaries - element = 0" << G4endl;
337  return;
338  }
339 
340 
341  G4IonisParamElm* ionisation = element->GetIonisation();
342 
343  if (ionisation == 0)
344  {
345  G4cout << "G4LivermorePolarizedGammaConversionModel::SampleSecondaries - ionisation = 0" << G4endl;
346  return;
347  }
348 
349  // Extract Coulomb factor for this Element
350 
351  G4double fZ = 8. * (ionisation->GetlogZ3());
352  if (photonEnergy > 50. * MeV) fZ += 8. * (element->GetfCoulomb());
353 
354  // Limits of the screening variable
355  G4double screenFactor = 136. * epsilon0Local / (element->GetIonisation()->GetZ3()) ;
356  G4double screenMax = G4Exp ((42.24 - fZ)/8.368) - 0.952 ;
357  G4double screenMin = std::min(4.*screenFactor,screenMax) ;
358 
359  // Limits of the energy sampling
360  G4double epsilon1 = 0.5 - 0.5 * sqrt(1. - screenMin / screenMax) ;
361  G4double epsilonMin = std::max(epsilon0Local,epsilon1);
362  G4double epsilonRange = 0.5 - epsilonMin ;
363 
364  // Sample the energy rate of the created electron (or positron)
365  G4double screen;
366  G4double gReject ;
367 
368  G4double f10 = ScreenFunction1(screenMin) - fZ;
369  G4double f20 = ScreenFunction2(screenMin) - fZ;
370  G4double normF1 = std::max(f10 * epsilonRange * epsilonRange,0.);
371  G4double normF2 = std::max(1.5 * f20,0.);
372 
373  do {
374  if (normF1 / (normF1 + normF2) > G4UniformRand() )
375  {
376  epsilon = 0.5 - epsilonRange * pow(G4UniformRand(), 0.3333) ;
377  screen = screenFactor / (epsilon * (1. - epsilon));
378  gReject = (ScreenFunction1(screen) - fZ) / f10 ;
379  }
380  else
381  {
382  epsilon = epsilonMin + epsilonRange * G4UniformRand();
383  screen = screenFactor / (epsilon * (1 - epsilon));
384  gReject = (ScreenFunction2(screen) - fZ) / f20 ;
385 
386 
387  }
388  } while ( gReject < G4UniformRand() );
389 
390  } // End of epsilon sampling
391 
392  // Fix charges randomly
393 
394  G4double electronTotEnergy;
395  G4double positronTotEnergy;
396 
397 
398  // if (G4int(2*G4UniformRand()))
399  if (G4UniformRand() > 0.5)
400  {
401  electronTotEnergy = (1. - epsilon) * photonEnergy;
402  positronTotEnergy = epsilon * photonEnergy;
403  }
404  else
405  {
406  positronTotEnergy = (1. - epsilon) * photonEnergy;
407  electronTotEnergy = epsilon * photonEnergy;
408  }
409 
410  // Scattered electron (positron) angles. ( Z - axis along the parent photon)
411  // Universal distribution suggested by L. Urban (Geant3 manual (1993) Phys211),
412  // derived from Tsai distribution (Rev. Mod. Phys. 49, 421 (1977)
413 
414 /*
415  G4double u;
416  const G4double a1 = 0.625;
417  G4double a2 = 3. * a1;
418 
419  if (0.25 > G4UniformRand())
420  {
421  u = - log(G4UniformRand() * G4UniformRand()) / a1 ;
422  }
423  else
424  {
425  u = - log(G4UniformRand() * G4UniformRand()) / a2 ;
426  }
427 */
428 
429  G4double Ene = electronTotEnergy/electron_mass_c2; // Normalized energy
430 
431  G4double cosTheta = 0.;
432  G4double sinTheta = 0.;
433 
434  SetTheta(&cosTheta,&sinTheta,Ene);
435 
436  // G4double theta = u * electron_mass_c2 / photonEnergy ;
437  // G4double phi = twopi * G4UniformRand() ;
438 
439  G4double phi,psi=0.;
440 
441  //corrected e+ e- angular angular distribution //preliminary!
442 
443  // if(photonEnergy>50*MeV)
444  // {
445  phi = SetPhi(photonEnergy);
446  psi = SetPsi(photonEnergy,phi);
447  // }
448  //else
449  // {
450  //psi = G4UniformRand()*2.*pi;
451  //phi = pi; // coplanar
452  // }
453 
454  Psi = psi;
455  Phi = phi;
456  //G4cout << "PHI " << phi << G4endl;
457  //G4cout << "PSI " << psi << G4endl;
458 
459  G4double phie, phip;
460  G4double choice, choice2;
461  choice = G4UniformRand();
462  choice2 = G4UniformRand();
463 
464  if (choice2 <= 0.5)
465  {
466  // do nothing
467  // phi = phi;
468  }
469  else
470  {
471  phi = -phi;
472  }
473 
474  if (choice <= 0.5)
475  {
476  phie = psi; //azimuthal angle for the electron
477  phip = phie+phi; //azimuthal angle for the positron
478  }
479  else
480  {
481  // opzione 1 phie / phip equivalenti
482 
483  phip = psi; //azimuthal angle for the positron
484  phie = phip + phi; //azimuthal angle for the electron
485  }
486 
487 
488  // Electron Kinematics
489 
490  G4double dirX = sinTheta*cos(phie);
491  G4double dirY = sinTheta*sin(phie);
492  G4double dirZ = cosTheta;
493  G4ThreeVector electronDirection(dirX,dirY,dirZ);
494 
495  // Kinematics of the created pair:
496  // the electron and positron are assumed to have a symetric angular
497  // distribution with respect to the Z axis along the parent photon
498 
499  //G4double localEnergyDeposit = 0. ;
500 
501  G4double electronKineEnergy = std::max(0.,electronTotEnergy - electron_mass_c2) ;
502 
503  SystemOfRefChange(gammaDirection0,electronDirection,
504  gammaPolarization0);
505 
507  electronDirection,
508  electronKineEnergy);
509 
510  // The e+ is always created (even with kinetic energy = 0) for further annihilation
511 
512  Ene = positronTotEnergy/electron_mass_c2; // Normalized energy
513 
514  cosTheta = 0.;
515  sinTheta = 0.;
516 
517  SetTheta(&cosTheta,&sinTheta,Ene);
518 
519  // Positron Kinematics
520 
521  dirX = sinTheta*cos(phip);
522  dirY = sinTheta*sin(phip);
523  dirZ = cosTheta;
524  G4ThreeVector positronDirection(dirX,dirY,dirZ);
525 
526  G4double positronKineEnergy = std::max(0.,positronTotEnergy - electron_mass_c2) ;
527  SystemOfRefChange(gammaDirection0,positronDirection,
528  gammaPolarization0);
529 
530  // Create G4DynamicParticle object for the particle2
532  positronDirection, positronKineEnergy);
533 
534 
535  fvect->push_back(particle1);
536  fvect->push_back(particle2);
537 
538  // Kill the incident photon
539  fParticleChange->SetProposedKineticEnergy(0.);
540  fParticleChange->ProposeTrackStatus(fStopAndKill);
541 
542 }
G4double GetKineticEnergy() const
G4double GetfCoulomb() const
Definition: G4Element.hh:191
G4ParticleDefinition * GetDefinition() const
bool isOrthogonal(const Hep3Vector &v, double epsilon=tolerance) const
Definition: SpaceVector.cc:237
double howOrthogonal(const Hep3Vector &v) const
Definition: SpaceVector.cc:219
static constexpr double electron_mass_c2
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
G4double GetlogZ3() const
G4double G4Exp(G4double initial_x)
Exponential Function double precision.
Definition: G4Exp.hh:183
static G4Positron * Positron()
Definition: G4Positron.cc:94
T max(const T t1, const T t2)
brief Return the largest of the two arguments
G4IonisParamElm * GetIonisation() const
Definition: G4Element.hh:199
const G4ThreeVector & GetPolarization() const
T min(const T t1, const T t2)
brief Return the smallest of the two arguments
static G4Electron * Electron()
Definition: G4Electron.cc:94
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
double mag() const
G4double GetZ3() const
double epsilon(double density, double temperature)
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:541

Here is the call graph for this function:


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