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

#include <G4LivermorePolarizedGammaConversionModel.hh>

+ Inheritance 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 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 42 of file G4LivermorePolarizedGammaConversionModel.hh.

Constructor & Destructor Documentation

◆ G4LivermorePolarizedGammaConversionModel()

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

Definition at line 50 of file G4LivermorePolarizedGammaConversionModel.cc.

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

◆ ~G4LivermorePolarizedGammaConversionModel()

G4LivermorePolarizedGammaConversionModel::~G4LivermorePolarizedGammaConversionModel ( )
virtual

Definition at line 75 of file G4LivermorePolarizedGammaConversionModel.cc.

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

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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 220 of file G4LivermorePolarizedGammaConversionModel.cc.

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

◆ Initialise()

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

Implements G4VEmModel.

Definition at line 89 of file G4LivermorePolarizedGammaConversionModel.cc.

91{
92 if (verboseLevel > 1)
93 {
94 G4cout << "Calling1 G4LivermorePolarizedGammaConversionModel::Initialise()"
95 << G4endl
96 << "Energy range: "
97 << LowEnergyLimit() / MeV << " MeV - "
98 << HighEnergyLimit() / GeV << " GeV"
99 << G4endl;
100 }
101
102
103 if(IsMaster())
104 {
105
106 // Initialise element selector
107
108 InitialiseElementSelectors(particle, cuts);
109
110 // Access to elements
111
112 char* path = std::getenv("G4LEDATA");
113
114 G4ProductionCutsTable* theCoupleTable =
116
117 G4int numOfCouples = theCoupleTable->GetTableSize();
118
119 for(G4int i=0; i<numOfCouples; ++i)
120 {
121 const G4Material* material =
122 theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
123 const G4ElementVector* theElementVector = material->GetElementVector();
124 G4int nelm = material->GetNumberOfElements();
125
126 for (G4int j=0; j<nelm; ++j)
127 {
128 G4int Z = (G4int)(*theElementVector)[j]->GetZ();
129 if(Z < 1) { Z = 1; }
130 else if(Z > maxZ) { Z = maxZ; }
131 if(!data[Z]) { ReadData(Z, path); }
132 }
133 }
134 }
135 if(isInitialised) { return; }
136 fParticleChange = GetParticleChangeForGamma();
137 isInitialised = true;
138
139}
std::vector< G4Element * > G4ElementVector
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()
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:133
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:652
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:645
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:148

◆ InitialiseForElement()

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

Reimplemented from G4VEmModel.

Definition at line 1086 of file G4LivermorePolarizedGammaConversionModel.cc.

1089{
1090 G4AutoLock l(&LivermorePolarizedGammaConversionModelMutex);
1091 // G4cout << "G4LivermorePolarizedGammaConversionModel::InitialiseForElement Z= "
1092 // << Z << G4endl;
1093 if(!data[Z]) { ReadData(Z); }
1094 l.unlock();
1095}

Referenced by ComputeCrossSectionPerAtom().

◆ InitialiseLocal()

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

Reimplemented from G4VEmModel.

Definition at line 144 of file G4LivermorePolarizedGammaConversionModel.cc.

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

◆ MinPrimaryEnergy()

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

Reimplemented from G4VEmModel.

Definition at line 152 of file G4LivermorePolarizedGammaConversionModel.cc.

154{
155 return lowEnergyLimit;
156}

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 268 of file G4LivermorePolarizedGammaConversionModel.cc.

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

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