Geant4  10.03.p03
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
G4BoldyshevTripletModel Class Reference

#include <G4BoldyshevTripletModel.hh>

Inheritance diagram for G4BoldyshevTripletModel:
Collaboration diagram for G4BoldyshevTripletModel:

Public Member Functions

 G4BoldyshevTripletModel (const G4ParticleDefinition *p=0, const G4String &nam="BoldyshevTripletConversion")
 
virtual ~G4BoldyshevTripletModel ()
 
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 SetFluctuationFlag (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
 
G4bool lossFlucFlag
 
- Static Protected Attributes inherited from G4VEmModel
static const G4double inveplus = 1.0/CLHEP::eplus
 

Detailed Description

Definition at line 38 of file G4BoldyshevTripletModel.hh.

Constructor & Destructor Documentation

G4BoldyshevTripletModel::G4BoldyshevTripletModel ( const G4ParticleDefinition p = 0,
const G4String nam = "BoldyshevTripletConversion" 
)

Definition at line 46 of file G4BoldyshevTripletModel.cc.

47  :G4VEmModel(nam),isInitialised(false),smallEnergy(4.*MeV)
48 {
49  fParticleChange = 0;
50 
51  lowEnergyLimit = 4.0*electron_mass_c2;
52 
53  verboseLevel= 0;
54  // Verbosity scale for debugging purposes:
55  // 0 = nothing
56  // 1 = calculation of cross sections, file openings...
57  // 2 = entering in methods
58 
59  if(verboseLevel > 0)
60  {
61  G4cout << "G4BoldyshevTripletModel is constructed " << G4endl;
62  }
63 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:68
G4GLOB_DLL std::ostream G4cout
float electron_mass_c2
Definition: hepunit.py:274
#define G4endl
Definition: G4ios.hh:61
static constexpr double MeV
Definition: G4SIunits.hh:214
G4BoldyshevTripletModel::~G4BoldyshevTripletModel ( )
virtual

Definition at line 67 of file G4BoldyshevTripletModel.cc.

68 {
69  if(IsMaster()) {
70  for(G4int i=0; i<maxZ; ++i) {
71  if(data[i]) {
72  delete data[i];
73  data[i] = 0;
74  }
75  }
76  }
77 }
G4bool IsMaster() const
Definition: G4VEmModel.hh:725
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 G4BoldyshevTripletModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  kinEnergy,
G4double  Z,
G4double  A = 0,
G4double  cut = 0,
G4double  emax = DBL_MAX 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 214 of file G4BoldyshevTripletModel.cc.

218 {
219  if (verboseLevel > 1)
220  {
221  G4cout << "Calling ComputeCrossSectionPerAtom() of G4BoldyshevTripletModel"
222  << G4endl;
223  }
224 
225  if (GammaEnergy < lowEnergyLimit) { return 0.0; }
226 
227  G4double xs = 0.0;
228 
229  G4int intZ=G4int(Z);
230 
231  if(intZ < 1 || intZ > maxZ) { return xs; }
232 
233  G4LPhysicsFreeVector* pv = data[intZ];
234 
235  // if element was not initialised
236  // do initialisation safely for MT mode
237  if(!pv)
238  {
239  InitialiseForElement(0, intZ);
240  pv = data[intZ];
241  if(!pv) { return xs; }
242  }
243  // x-section is taken from the table
244  xs = pv->Value(GammaEnergy);
245 
246  if(verboseLevel > 0)
247  {
248  G4int n = pv->GetVectorLength() - 1;
249  G4cout << "****** DEBUG: tcs value for Z=" << Z << " at energy (MeV)="
250  << GammaEnergy/MeV << G4endl;
251  G4cout << " cs (Geant4 internal unit)=" << xs << G4endl;
252  G4cout << " -> first cs value in EADL data file (iu) =" << (*pv)[0] << G4endl;
253  G4cout << " -> last cs value in EADL data file (iu) =" << (*pv)[n] << G4endl;
254  G4cout << "*********************************************************" << G4endl;
255  }
256 
257  return xs;
258 
259 }
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
const G4int n
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
virtual void InitialiseForElement(const G4ParticleDefinition *, G4int Z)

Here is the call graph for this function:

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

Implements G4VEmModel.

Definition at line 81 of file G4BoldyshevTripletModel.cc.

84 {
85  if (verboseLevel > 1)
86  {
87  G4cout << "Calling Initialise() of G4BoldyshevTripletModel."
88  << G4endl
89  << "Energy range: "
90  << LowEnergyLimit() / MeV << " MeV - "
91  << HighEnergyLimit() / GeV << " GeV"
92  << G4endl;
93  }
94 
95  if(IsMaster())
96  {
97 
98  // Initialise element selector
99 
100  InitialiseElementSelectors(particle, cuts);
101 
102  // Access to elements
103 
104  char* path = getenv("G4LEDATA");
105 
106  G4ProductionCutsTable* theCoupleTable =
108 
109  G4int numOfCouples = theCoupleTable->GetTableSize();
110 
111  for(G4int i=0; i<numOfCouples; ++i)
112  {
113  const G4Material* material =
114  theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
115  const G4ElementVector* theElementVector = material->GetElementVector();
116  G4int nelm = material->GetNumberOfElements();
117 
118  for (G4int j=0; j<nelm; ++j)
119  {
120  G4int Z = (G4int)(*theElementVector)[j]->GetZ();
121  if(Z < 1) { Z = 1; }
122  else if(Z > maxZ) { Z = maxZ; }
123  if(!data[Z]) { ReadData(Z, path); }
124  }
125  }
126  }
127  if(isInitialised) { return; }
128  fParticleChange = GetParticleChangeForGamma();
129  isInitialised = true;
130 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:643
std::vector< G4Element * > G4ElementVector
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:146
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:636
G4bool IsMaster() const
Definition: G4VEmModel.hh:725
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
string material
Definition: eplot.py:19
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 G4BoldyshevTripletModel::InitialiseForElement ( const G4ParticleDefinition ,
G4int  Z 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 478 of file G4BoldyshevTripletModel.cc.

481 {
482  G4AutoLock l(&BoldyshevTripletModelMutex);
483  // G4cout << "G4BoldyshevTripletModel::InitialiseForElement Z= "
484  // << Z << G4endl;
485  if(!data[Z]) { ReadData(Z); }
486  l.unlock();
487 }
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 G4BoldyshevTripletModel::InitialiseLocal ( const G4ParticleDefinition ,
G4VEmModel masterModel 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 134 of file G4BoldyshevTripletModel.cc.

136 {
138 }
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:809
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:817

Here is the call graph for this function:

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

Reimplemented from G4VEmModel.

Definition at line 143 of file G4BoldyshevTripletModel.cc.

146 {
147  return lowEnergyLimit;
148 }
void G4BoldyshevTripletModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple ,
const G4DynamicParticle aDynamicGamma,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 263 of file G4BoldyshevTripletModel.cc.

268 {
269 
270  // The energies of the secondary particles are sampled using // a modified Wheeler-Lamb model (see PhysRevD 7 (1973), 26)
271 
272  if (verboseLevel > 1) {
273  G4cout << "Calling SampleSecondaries() of G4BoldyshevTripletModel"
274  << G4endl;
275  }
276 
277  G4double photonEnergy = aDynamicGamma->GetKineticEnergy();
278  G4ParticleMomentum photonDirection = aDynamicGamma->GetMomentumDirection();
279 
280  G4double epsilon ;
281  // G4double epsilon0Local = electron_mass_c2 / photonEnergy ;
282 
283 
285  G4double positronTotEnergy, electronTotEnergy, thetaEle, thetaPos;
286  G4double ener_re=0., theta_re, phi_re, phi;
287 
288  // Calculo de theta - elecron de recoil
289 
290  G4double energyThreshold = sqrt(2.)*electron_mass_c2; // -> momentumThreshold_N = 1
291 
292  energyThreshold = 1.1*electron_mass_c2;
293  // G4cout << energyThreshold << G4endl;
294 
295 
296  G4double momentumThreshold_c = sqrt(energyThreshold * energyThreshold - electron_mass_c2*electron_mass_c2); // momentun in MeV/c unit
297  G4double momentumThreshold_N = momentumThreshold_c/electron_mass_c2; // momentun in mc unit
298 
299  // Calculation of recoil electron production
300 
301 
302  G4double SigmaTot = (28./9.) * std::log ( 2.* photonEnergy / electron_mass_c2 ) - 218. / 27. ;
303  G4double X_0 = 2. * ( sqrt(momentumThreshold_N*momentumThreshold_N + 1) -1 );
304  G4double SigmaQ = (82./27. - (14./9.) * log (X_0) + 4./15.*X_0 - 0.0348 * X_0 * X_0);
305  G4double recoilProb = G4UniformRand();
306  //G4cout << "SIGMA TOT " << SigmaTot << " " << "SigmaQ " << SigmaQ << " " << SigmaQ/SigmaTot << " " << recoilProb << G4endl;
307 
308 
309  if (recoilProb >= SigmaQ/SigmaTot) // create electron recoil
310  {
311 
312  G4double cosThetaMax = ( ( energyThreshold - electron_mass_c2 ) / (momentumThreshold_c) + electron_mass_c2*
313  ( energyThreshold + electron_mass_c2 ) / (photonEnergy*momentumThreshold_c) );
314 
315 
316  if (cosThetaMax > 1) G4cout << "ERRORE " << G4endl;
317 
318  G4double r1;
319  G4double r2;
320  G4double are, bre, loga, f1_re, greject, cost;
321 
322  do {
323  r1 = G4UniformRand();
324  r2 = G4UniformRand();
325  // cost = (pow(4./enern,0.5*r1)) ;
326 
327  cost = pow(cosThetaMax,r1);
328  theta_re = acos(cost);
329  are = 1./(14.*cost*cost);
330  bre = (1.-5.*cost*cost)/(2.*cost);
331  loga = log((1.+ cost)/(1.- cost));
332  f1_re = 1. - bre*loga;
333 
334  if ( theta_re >= 4.47*CLHEP::pi/180.)
335  {
336  greject = are*f1_re;
337  } else {
338  greject = 1. ;
339  }
340  } while(greject < r2);
341 
342  // Calculo de phi - elecron de recoil
343 
344  G4double r3, r4, rt;
345 
346  do {
347  r3 = G4UniformRand();
348  r4 = G4UniformRand();
349  phi_re = twopi*r3 ;
350  G4double sint2 = 1. - cost*cost ;
351  G4double fp = 1. - sint2*loga/(2.*cost) ;
352  rt = (1.-cos(2.*phi_re)*fp/f1_re)/(2.*pi) ;
353 
354  } while(rt < r4);
355 
356  // Calculo de la energia - elecron de recoil - relacion momento maximo <-> angulo
357 
358  G4double S = electron_mass_c2*(2.* photonEnergy + electron_mass_c2);
359 
362  *(S - electron_mass_c2*electron_mass_c2)*sin(theta_re)*sin(theta_re);
363  ener_re = electron_mass_c2 * (S + electron_mass_c2*electron_mass_c2)/sqrt(D2);
364 
365  // New Recoil energy calculation
366 
367  G4double momentum_recoil = 2* (electron_mass_c2) * (std::cos(theta_re)/(std::sin(phi_re)*std::sin(phi_re)));
368  G4double ener_recoil = sqrt( momentum_recoil*momentum_recoil + electron_mass_c2*electron_mass_c2);
369  ener_re = ener_recoil;
370 
371  // G4cout << "electron de retroceso " << ener_re << " " << theta_re << " " << phi_re << G4endl;
372 
373  // Recoil electron creation
374  G4double dxEle_re=sin(theta_re)*std::cos(phi_re),dyEle_re=sin(theta_re)*std::sin(phi_re), dzEle_re=cos(theta_re);
375 
376  G4double electronRKineEnergy = std::max(0.,ener_re - electron_mass_c2) ;
377 
378  G4ThreeVector electronRDirection (dxEle_re, dyEle_re, dzEle_re);
379  electronRDirection.rotateUz(photonDirection);
380 
382  electronRDirection,
383  electronRKineEnergy);
384  fvect->push_back(particle3);
385 
386  }
387  else
388  {
389  // deposito la energia ener_re - electron_mass_c2
390  // G4cout << "electron de retroceso " << ener_re << G4endl;
391 
392  fParticleChange->ProposeLocalEnergyDeposit(ener_re - electron_mass_c2);
393  }
394 
395 
396  // Depaola (2004) suggested distribution for e+e- energy
397 
398  // G4double t = 0.5*asinh(momentumThreshold_N);
399  G4double t = 0.5*log(momentumThreshold_N + sqrt(momentumThreshold_N*momentumThreshold_N+1));
400  //G4cout << 0.5*asinh(momentumThreshold_N) << " " << t << G4endl;
401  G4double J1 = 0.5*(t*cosh(t)/sinh(t) - log(2.*sinh(t)));
402  G4double J2 = (-2./3.)*log(2.*sinh(t)) + t*cosh(t)/sinh(t) + (sinh(t)-t*pow(cosh(t),3))/(3.*pow(sinh(t),3));
403  G4double b = 2.*(J1-J2)/J1;
404 
405  G4double n = 1 - b/6.;
406  G4double re=0.;
407  re = G4UniformRand();
408  G4double a = 0.;
409 
410  G4double b1 = 16. - 3.*b - 36.*b*re*n + 36.*b*pow(re,2.)*pow(n,2.) +
411  6.*pow(b,2.)*re*n;
412  a = pow((b1/b),0.5);
413  G4double c1 = (-6. + 12.*re*n + b + 2*a)*pow(b,2.);
414  epsilon = (pow(c1,1./3.))/(2.*b) + (b-4.)/(2.*pow(c1,1./3.))+0.5;
415 
416  G4double photonEnergy1 = photonEnergy - ener_re ; // resto al foton la energia del electron de retro.
417  positronTotEnergy = epsilon*photonEnergy1;
418  electronTotEnergy = photonEnergy1 - positronTotEnergy; // temporarly
419 
420  G4double momento_e = sqrt(electronTotEnergy*electronTotEnergy -
422  G4double momento_p = sqrt(positronTotEnergy*positronTotEnergy -
424 
425  thetaEle = acos((sqrt(p0*p0/(momento_e*momento_e) +1.)- p0/momento_e)) ;
426  thetaPos = acos((sqrt(p0*p0/(momento_p*momento_p) +1.)- p0/momento_p)) ;
427  phi = twopi * G4UniformRand();
428 
429  G4double dxEle= std::sin(thetaEle)*std::cos(phi),dyEle= std::sin(thetaEle)*std::sin(phi),dzEle=std::cos(thetaEle);
430  G4double dxPos=-std::sin(thetaPos)*std::cos(phi),dyPos=-std::sin(thetaPos)*std::sin(phi),dzPos=std::cos(thetaPos);
431 
432  // Kinematics of the created pair:
433 
434  // the electron and positron are assumed to have a symetric angular
435  // distribution with respect to the Z axis along the parent photon
436 
437 
438  G4double electronKineEnergy = std::max(0.,electronTotEnergy - electron_mass_c2) ;
439 
440 
441  G4ThreeVector electronDirection (dxEle, dyEle, dzEle);
442  electronDirection.rotateUz(photonDirection);
443 
445  electronDirection,
446  electronKineEnergy);
447 
448  // The e+ is always created (even with kinetic energy = 0) for further annihilation
449 
450  G4double positronKineEnergy = std::max(0.,positronTotEnergy - electron_mass_c2) ;
451 
452  G4ThreeVector positronDirection (dxPos, dyPos, dzPos);
453  positronDirection.rotateUz(photonDirection);
454 
455  // Create G4DynamicParticle object for the particle2
456 
458  positronDirection, positronKineEnergy);
459  // Fill output vector
460 
461  fvect->push_back(particle1);
462  fvect->push_back(particle2);
463 
464 
465  // kill incident photon
466  fParticleChange->SetProposedKineticEnergy(0.);
467  fParticleChange->ProposeTrackStatus(fStopAndKill);
468 
469 
470 
471 }
double S(double temp)
G4double GetKineticEnergy() const
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
void ProposeLocalEnergyDeposit(G4double anEnergyPart)
static constexpr double twopi
Definition: G4SIunits.hh:76
tuple b
Definition: test.py:12
#define G4UniformRand()
Definition: Randomize.hh:97
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
float electron_mass_c2
Definition: hepunit.py:274
const G4int n
static G4Positron * Positron()
Definition: G4Positron.cc:94
T max(const T t1, const T t2)
brief Return the largest of the two arguments
static G4Electron * Electron()
Definition: G4Electron.cc:94
void SetProposedKineticEnergy(G4double proposedKinEnergy)
#define G4endl
Definition: G4ios.hh:61
static constexpr double pi
Definition: G4SIunits.hh:75
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
double epsilon(double density, double temperature)
static constexpr double pi
Definition: SystemOfUnits.h:54
tuple c1
Definition: plottest35.py:14

Here is the call graph for this function:


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