Geant4 11.1.1
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=nullptr, const G4String &nam="LivermorePolarizedGammaConversion")
 
virtual ~G4LivermorePolarizedGammaConversionModel ()
 
void Initialise (const G4ParticleDefinition *, const G4DataVector &) override
 
void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel) override
 
void InitialiseForElement (const G4ParticleDefinition *, G4int Z) override
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX) override
 
void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy) override
 
G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double) override
 
G4LivermorePolarizedGammaConversionModeloperator= (const G4LivermorePolarizedGammaConversionModel &right)=delete
 
 G4LivermorePolarizedGammaConversionModel (const G4LivermorePolarizedGammaConversionModel &)=delete
 
- 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 *, const G4double &length, G4double &eloss)
 
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 FillNumberOfSecondaries (G4int &numberOfTriplets, G4int &numberOfRecoil)
 
virtual void ModelDescription (std::ostream &outFile) const
 
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)
 
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)
 
const G4ElementGetCurrentElement (const G4Material *mat=nullptr) const
 
G4int SelectRandomAtomNumber (const G4Material *) const
 
const G4IsotopeGetCurrentIsotope (const G4Element *elm=nullptr) const
 
G4int SelectIsotopeNumber (const G4Element *) const
 
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 *)
 
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 = nullptr
 
G4VParticleChangepParticleChange = nullptr
 
G4PhysicsTablexSectionTable = nullptr
 
const G4MaterialpBaseMaterial = nullptr
 
const std::vector< G4double > * theDensityFactor = nullptr
 
const std::vector< G4int > * theDensityIdx = nullptr
 
G4double inveplus
 
G4double pFactor = 1.0
 
size_t currentCoupleIndex = 0
 
size_t basedCoupleIndex = 0
 
G4bool lossFlucFlag = true
 

Detailed Description

Definition at line 38 of file G4LivermorePolarizedGammaConversionModel.hh.

Constructor & Destructor Documentation

◆ G4LivermorePolarizedGammaConversionModel() [1/2]

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

Definition at line 51 of file G4LivermorePolarizedGammaConversionModel.cc.

53 :G4VEmModel(nam), smallEnergy(2.*MeV), isInitialised(false)
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}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout

◆ ~G4LivermorePolarizedGammaConversionModel()

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] = nullptr;
83 }
84 }
85 }
86}
int G4int
Definition: G4Types.hh:85
G4bool IsMaster() const
Definition: G4VEmModel.hh:725

◆ G4LivermorePolarizedGammaConversionModel() [2/2]

G4LivermorePolarizedGammaConversionModel::G4LivermorePolarizedGammaConversionModel ( const G4LivermorePolarizedGammaConversionModel )
delete

Member Function Documentation

◆ ComputeCrossSectionPerAtom()

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

Reimplemented from G4VEmModel.

Definition at line 211 of file G4LivermorePolarizedGammaConversionModel.cc.

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

◆ Initialise()

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

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 if(IsMaster())
104 {
105 // Initialise element selector
106 InitialiseElementSelectors(particle, cuts);
107
108 // Access to elements
109 const char* path = G4FindDataDir("G4LEDATA");
110
111 G4ProductionCutsTable* theCoupleTable =
113
114 G4int numOfCouples = (G4int)theCoupleTable->GetTableSize();
115
116 for(G4int i=0; i<numOfCouples; ++i)
117 {
118 const G4Material* material =
119 theCoupleTable->GetMaterialCutsCouple(i)->GetMaterial();
120 const G4ElementVector* theElementVector = material->GetElementVector();
121 std::size_t nelm = material->GetNumberOfElements();
122
123 for (std::size_t j=0; j<nelm; ++j)
124 {
125 G4int Z = (G4int)(*theElementVector)[j]->GetZ();
126 if(Z < 1) { Z = 1; }
127 else if(Z > maxZ) { Z = maxZ; }
128 if(!data[Z]) { ReadData(Z, path); }
129 }
130 }
131 }
132 if(isInitialised) { return; }
133 fParticleChange = GetParticleChangeForGamma();
134 isInitialised = true;
135}
std::vector< const G4Element * > G4ElementVector
const char * G4FindDataDir(const char *)
const G4Material * GetMaterial() const
const G4ElementVector * GetElementVector() const
Definition: G4Material.hh:185
size_t GetNumberOfElements() const
Definition: G4Material.hh:181
const G4MaterialCutsCouple * GetMaterialCutsCouple(G4int i) const
std::size_t GetTableSize() const
static G4ProductionCutsTable * GetProductionCutsTable()
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:124
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:641
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:634
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:139

◆ InitialiseForElement()

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

Reimplemented from G4VEmModel.

Definition at line 984 of file G4LivermorePolarizedGammaConversionModel.cc.

987{
988 G4AutoLock l(&LivermorePolarizedGammaConversionModelMutex);
989 if(!data[Z]) { ReadData(Z); }
990 l.unlock();
991}

Referenced by ComputeCrossSectionPerAtom().

◆ InitialiseLocal()

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

Reimplemented from G4VEmModel.

Definition at line 139 of file G4LivermorePolarizedGammaConversionModel.cc.

141{
143}
void SetElementSelectors(std::vector< G4EmElementSelector * > *)
Definition: G4VEmModel.hh:831
std::vector< G4EmElementSelector * > * GetElementSelectors()
Definition: G4VEmModel.hh:823

◆ MinPrimaryEnergy()

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

Reimplemented from G4VEmModel.

Definition at line 147 of file G4LivermorePolarizedGammaConversionModel.cc.

149{
150 return lowEnergyLimit;
151}

◆ operator=()

G4LivermorePolarizedGammaConversionModel & G4LivermorePolarizedGammaConversionModel::operator= ( const G4LivermorePolarizedGammaConversionModel right)
delete

◆ SampleSecondaries()

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

Implements G4VEmModel.

Definition at line 259 of file G4LivermorePolarizedGammaConversionModel.cc.

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

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