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

#include <G4LivermoreGammaConversionModelRC.hh>

Inheritance diagram for G4LivermoreGammaConversionModelRC:
G4VEmModel

Public Member Functions

 G4LivermoreGammaConversionModelRC (const G4ParticleDefinition *p=0, const G4String &nam="LivermoreConversion")
 
virtual ~G4LivermoreGammaConversionModelRC ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A=0, G4double cut=0, G4double emax=DBL_MAX)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
virtual void 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 ChargeSquareRatio (const G4Track &)
 
virtual G4double GetChargeSquareRatio (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual G4double GetParticleCharge (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void StartTracking (G4Track *)
 
virtual void CorrectionsAlongStep (const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double &eloss, G4double &niel, G4double length)
 
virtual G4double Value (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy)
 
virtual G4double MinPrimaryEnergy (const G4Material *, const G4ParticleDefinition *, G4double cut=0.0)
 
virtual G4double MinEnergyCut (const G4ParticleDefinition *, const G4MaterialCutsCouple *)
 
virtual void SetupForMaterial (const G4ParticleDefinition *, const G4Material *, G4double kineticEnergy)
 
virtual void DefineForRegion (const G4Region *)
 
void InitialiseElementSelectors (const G4ParticleDefinition *, const G4DataVector &)
 
std::vector
< G4EmElementSelector * > * 
GetElementSelectors ()
 
void SetElementSelectors (std::vector< G4EmElementSelector * > *)
 
G4double ComputeDEDX (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=DBL_MAX)
 
G4double CrossSection (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeMeanFreePath (const G4ParticleDefinition *, G4double kineticEnergy, const G4Material *, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, const G4Element *, G4double kinEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectIsotopeNumber (const G4Element *)
 
const G4ElementSelectRandomAtom (const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
const G4ElementSelectRandomAtom (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
 
G4int SelectRandomAtomNumber (const G4Material *)
 
void SetParticleChange (G4VParticleChange *, G4VEmFluctuationModel *f=0)
 
void SetCrossSectionTable (G4PhysicsTable *, G4bool isLocal)
 
G4ElementDataGetElementData ()
 
G4PhysicsTableGetCrossSectionTable ()
 
G4VEmFluctuationModelGetModelOfFluctuations ()
 
G4VEmAngularDistributionGetAngularDistribution ()
 
void SetAngularDistribution (G4VEmAngularDistribution *)
 
G4double HighEnergyLimit () const
 
G4double LowEnergyLimit () const
 
G4double HighEnergyActivationLimit () const
 
G4double LowEnergyActivationLimit () const
 
G4double PolarAngleLimit () const
 
G4double SecondaryThreshold () const
 
G4bool LPMFlag () const
 
G4bool DeexcitationFlag () const
 
G4bool ForceBuildTableFlag () const
 
G4bool UseAngularGeneratorFlag () const
 
void SetAngularGeneratorFlag (G4bool)
 
void SetHighEnergyLimit (G4double)
 
void SetLowEnergyLimit (G4double)
 
void SetActivationHighEnergyLimit (G4double)
 
void SetActivationLowEnergyLimit (G4double)
 
G4bool IsActive (G4double kinEnergy)
 
void SetPolarAngleLimit (G4double)
 
void SetSecondaryThreshold (G4double)
 
void SetLPMFlag (G4bool val)
 
void SetDeexcitationFlag (G4bool val)
 
void SetForceBuildTable (G4bool val)
 
void SetMasterThread (G4bool val)
 
G4bool IsMaster () const
 
G4double MaxSecondaryKinEnergy (const G4DynamicParticle *dynParticle)
 
const G4StringGetName () const
 
void SetCurrentCouple (const G4MaterialCutsCouple *)
 
const G4ElementGetCurrentElement () const
 

Protected Member Functions

G4double GetMeanFreePath (const G4Track &aTrack, G4double previousStepSize, G4ForceCondition *condition)
 
- Protected Member Functions inherited from G4VEmModel
G4ParticleChangeForLossGetParticleChangeForLoss ()
 
G4ParticleChangeForGammaGetParticleChangeForGamma ()
 
virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kineticEnergy)
 
const G4MaterialCutsCoupleCurrentCouple () const
 
void SetCurrentElement (const G4Element *)
 

Protected Attributes

G4ParticleChangeForGammafParticleChange
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

Detailed Description

Definition at line 44 of file G4LivermoreGammaConversionModelRC.hh.

Constructor & Destructor Documentation

G4LivermoreGammaConversionModelRC::G4LivermoreGammaConversionModelRC ( const G4ParticleDefinition p = 0,
const G4String nam = "LivermoreConversion" 
)

Definition at line 50 of file G4LivermoreGammaConversionModelRC.cc.

References python.hepunit::electron_mass_c2, G4cout, G4endl, python.hepunit::GeV, python.hepunit::MeV, and G4VEmModel::SetHighEnergyLimit().

52  :G4VEmModel(nam),fParticleChange(0),smallEnergy(2.*MeV),isInitialised(false),
53  crossSectionHandler(0),meanFreePathTable(0)
54 {
55  lowEnergyLimit = 4.0*electron_mass_c2;
56  highEnergyLimit = 100 * GeV;
57  SetHighEnergyLimit(highEnergyLimit);
58 
59  verboseLevel= 0;
60  // Verbosity scale:
61  // 0 = nothing
62  // 1 = warning for energy non-conservation
63  // 2 = details of energy budget
64  // 3 = calculation of cross sections, file openings, sampling of atoms
65  // 4 = entering in methods
66 
67  if(verboseLevel > 0) {
68  G4cout << "Livermore Gamma conversion is constructed " << G4endl
69  << "Energy range: "
70  << lowEnergyLimit / MeV << " MeV - "
71  << highEnergyLimit / GeV << " GeV"
72  << G4endl;
73  }
74 }
G4VEmModel(const G4String &nam)
Definition: G4VEmModel.cc:65
void SetHighEnergyLimit(G4double)
Definition: G4VEmModel.hh:683
G4GLOB_DLL std::ostream G4cout
float electron_mass_c2
Definition: hepunit.py:274
#define G4endl
Definition: G4ios.hh:61
G4LivermoreGammaConversionModelRC::~G4LivermoreGammaConversionModelRC ( )
virtual

Definition at line 78 of file G4LivermoreGammaConversionModelRC.cc.

79 {
80  if (crossSectionHandler) delete crossSectionHandler;
81 }

Member Function Documentation

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

Reimplemented from G4VEmModel.

Definition at line 126 of file G4LivermoreGammaConversionModelRC.cc.

References G4VCrossSectionHandler::FindValue(), G4cout, and G4endl.

130 {
131  if (verboseLevel > 3) {
132  G4cout << "Calling ComputeCrossSectionPerAtom() of G4LivermoreGammaConversionModelRC"
133  << G4endl;
134  }
135  if (GammaEnergy < lowEnergyLimit || GammaEnergy > highEnergyLimit) return 0;
136 
137  G4double cs = crossSectionHandler->FindValue(G4int(Z), GammaEnergy);
138  return cs;
139 }
int G4int
Definition: G4Types.hh:78
G4double FindValue(G4int Z, G4double e) const
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76
G4double G4LivermoreGammaConversionModelRC::GetMeanFreePath ( const G4Track aTrack,
G4double  previousStepSize,
G4ForceCondition condition 
)
protected
void G4LivermoreGammaConversionModelRC::Initialise ( const G4ParticleDefinition ,
const G4DataVector  
)
virtual

Implements G4VEmModel.

Definition at line 86 of file G4LivermoreGammaConversionModelRC.cc.

References G4VCrossSectionHandler::Clear(), fParticleChange, G4cout, G4endl, G4VEmModel::GetParticleChangeForGamma(), python.hepunit::GeV, G4VEmModel::HighEnergyLimit(), G4VCrossSectionHandler::Initialise(), G4VCrossSectionHandler::LoadData(), G4VEmModel::LowEnergyLimit(), and python.hepunit::MeV.

88 {
89  if (verboseLevel > 3)
90  G4cout << "Calling G4LivermoreGammaConversionModelRC::Initialise()" << G4endl;
91 
92  if (crossSectionHandler)
93  {
94  crossSectionHandler->Clear();
95  delete crossSectionHandler;
96  }
97 
98  // Read data tables for all materials
99 
100  crossSectionHandler = new G4CrossSectionHandler();
101  crossSectionHandler->Initialise(0,lowEnergyLimit,100.*GeV,400);
102  G4String crossSectionFile = "pair/pp-cs-";
103  crossSectionHandler->LoadData(crossSectionFile);
104 
105  //
106 
107  if (verboseLevel > 2)
108  G4cout << "Loaded cross section files for Livermore Gamma Conversion model RC" << G4endl;
109 
110  if (verboseLevel > 0) {
111  G4cout << "Livermore Gamma Conversion model is initialized " << G4endl
112  << "Energy range: "
113  << LowEnergyLimit() / MeV << " MeV - "
114  << HighEnergyLimit() / GeV << " GeV"
115  << G4endl;
116  }
117 
118  if(isInitialised) return;
120  isInitialised = true;
121 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:599
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
G4GLOB_DLL std::ostream G4cout
void Initialise(G4VDataSetAlgorithm *interpolation=0, G4double minE=250 *CLHEP::eV, G4double maxE=100 *CLHEP::GeV, G4int numberOfBins=200, G4double unitE=CLHEP::MeV, G4double unitData=CLHEP::barn, G4int minZ=1, G4int maxZ=99)
void LoadData(const G4String &dataFile)
#define G4endl
Definition: G4ios.hh:61
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:121
void G4LivermoreGammaConversionModelRC::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle aDynamicGamma,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 143 of file G4LivermoreGammaConversionModelRC.cc.

References test::a, test::b, test::c, G4Electron::Electron(), python.hepunit::electron_mass_c2, fParticleChange, fStopAndKill, G4cout, G4endl, G4UniformRand, G4Gamma::Gamma(), G4DynamicParticle::GetDefinition(), G4Element::GetfCoulomb(), G4Element::GetIonisation(), G4DynamicParticle::GetKineticEnergy(), G4IonisParamElm::GetlogZ3(), G4DynamicParticle::GetMomentumDirection(), G4IonisParamElm::GetZ3(), G4INCL::Math::max(), python.hepunit::MeV, G4INCL::Math::min(), G4Positron::Positron(), G4VParticleChange::ProposeTrackStatus(), CLHEP::Hep3Vector::rotateUz(), G4VEmModel::SelectRandomAtom(), G4ParticleChangeForGamma::SetProposedKineticEnergy(), and python.hepunit::twopi.

148 {
149 
150 // The energies of the e+ e- secondaries are sampled using the Bethe - Heitler
151 // cross sections with Coulomb correction. A modified version of the random
152 // number techniques of Butcher & Messel is used (Nuc Phys 20(1960),15).
153 
154 // Note 1 : Effects due to the breakdown of the Born approximation at low
155 // energy are ignored.
156 // Note 2 : The differential cross section implicitly takes account of
157 // pair creation in both nuclear and atomic electron fields. However triplet
158 // prodution is not generated.
159 
160  if (verboseLevel > 3)
161  G4cout << "Calling SampleSecondaries() of G4LivermoreGammaConversionModelRC" << G4endl;
162 
163  G4double photonEnergy = aDynamicGamma->GetKineticEnergy();
164  G4ParticleMomentum photonDirection = aDynamicGamma->GetMomentumDirection();
165 
166  G4double epsilon ;
167  G4double epsilon0Local = electron_mass_c2 / photonEnergy ;
168  G4double electronTotEnergy;
169  G4double positronTotEnergy;
170 
171  G4double HardPhotonEnergy = 0.0;
172 
173 
174  // Do it fast if photon energy < 2. MeV
175  if (photonEnergy < smallEnergy )
176  {
177  epsilon = epsilon0Local + (0.5 - epsilon0Local) * G4UniformRand();
178 
179  if (G4int(2*G4UniformRand()))
180  {
181  electronTotEnergy = (1. - epsilon) * photonEnergy;
182  positronTotEnergy = epsilon * photonEnergy;
183  }
184  else
185  {
186  positronTotEnergy = (1. - epsilon) * photonEnergy;
187  electronTotEnergy = epsilon * photonEnergy;
188  }
189  }
190  else
191  {
192  // Select randomly one element in the current material
193  //const G4Element* element = crossSectionHandler->SelectRandomElement(couple,photonEnergy);
194  const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
195  const G4Element* element = SelectRandomAtom(couple,particle,photonEnergy);
196  G4cout << "G4LivermoreGammaConversionModelRC::SampleSecondaries" << G4endl;
197 
198  if (element == 0)
199  {
200  G4cout << "G4LivermoreGammaConversionModelRC::SampleSecondaries - element = 0"
201  << G4endl;
202  return;
203  }
204  G4IonisParamElm* ionisation = element->GetIonisation();
205  if (ionisation == 0)
206  {
207  G4cout << "G4LivermoreGammaConversionModelRC::SampleSecondaries - ionisation = 0"
208  << G4endl;
209  return;
210  }
211 
212  // Extract Coulomb factor for this Element
213  G4double fZ = 8. * (ionisation->GetlogZ3());
214  if (photonEnergy > 50. * MeV) fZ += 8. * (element->GetfCoulomb());
215 
216  // Limits of the screening variable
217  G4double screenFactor = 136. * epsilon0Local / (element->GetIonisation()->GetZ3()) ;
218  G4double screenMax = std::exp ((42.24 - fZ)/8.368) - 0.952 ;
219  G4double screenMin = std::min(4.*screenFactor,screenMax) ;
220 
221  // Limits of the energy sampling
222  G4double epsilon1 = 0.5 - 0.5 * std::sqrt(1. - screenMin / screenMax) ;
223  G4double epsilonMin = std::max(epsilon0Local,epsilon1);
224  G4double epsilonRange = 0.5 - epsilonMin ;
225 
226  // Sample the energy rate of the created electron (or positron)
227  G4double screen;
228  G4double gReject ;
229 
230  G4double f10 = ScreenFunction1(screenMin) - fZ;
231  G4double f20 = ScreenFunction2(screenMin) - fZ;
232  G4double normF1 = std::max(f10 * epsilonRange * epsilonRange,0.);
233  G4double normF2 = std::max(1.5 * f20,0.);
234  G4double a=393.3750918, b=115.3070201, c=810.6428451, d=19.96497475, e=1016.874592, f=1.936685510,
235  gLocal=751.2140962, h=0.099751048, i=299.9466339, j=0.002057250, k=49.81034926;
236  G4double aa=-18.6371131, bb=-1729.95248, cc=9450.971186, dd=106336.0145, ee=55143.09287, ff=-117602.840,
237  gg=-721455.467, hh=693957.8635, ii=156266.1085, jj=533209.9347;
238  G4double Rechazo = 0.;
239  G4double logepsMin = log(epsilonMin);
240  G4double NormaRC = a + b*logepsMin + c/logepsMin + d*pow(logepsMin,2.) + e/pow(logepsMin,2.) + f*pow(logepsMin,3.) +
241  gLocal/pow(logepsMin,3.) + h*pow(logepsMin,4.) + i/pow(logepsMin,4.) + j*pow(logepsMin,5.) +
242  k/pow(logepsMin,5.);
243 
244 
245  G4double HardPhotonThreshold = 0.08;
246  G4double r1, r2, r3, beta=0, gbeta, sigt = 582.068, sigh, rejet;
247  // , Pi = 2.*acos(0.);
248  G4double cg = (11./2.)/(exp(-11.*HardPhotonThreshold/2.)-exp(-11./2.));
249 
250 
251  r1 = G4UniformRand();
252  sigh = 1028.58*exp(-HardPhotonThreshold/0.09033) + 136.63; // sigma hard
253  if (r1 > 1.- sigh/sigt) {
254  r2 = G4UniformRand();
255  rejet = 0.;
256  while (r2 > rejet) {
257  r3 = G4UniformRand();
258  beta = (-2./11.)*log(exp(-0.08*11./2.)-r3*11./(2.*cg));
259  gbeta = exp(-11.*beta/2.);
260  rejet = fbeta(beta)/(8000.*gbeta);
261  }
262  HardPhotonEnergy = beta * photonEnergy;
263  }
264  else{
265  HardPhotonEnergy = 0.;
266  }
267 
268  photonEnergy -= HardPhotonEnergy;
269 
270  do {
271  do {
272  if (normF1 / (normF1 + normF2) > G4UniformRand() )
273  {
274  epsilon = 0.5 - epsilonRange * std::pow(G4UniformRand(), 0.3333) ;
275  screen = screenFactor / (epsilon * (1. - epsilon));
276  gReject = (ScreenFunction1(screen) - fZ) / f10 ;
277  }
278  else
279  {
280  epsilon = epsilonMin + epsilonRange * G4UniformRand();
281  screen = screenFactor / (epsilon * (1 - epsilon));
282  gReject = (ScreenFunction2(screen) - fZ) / f20 ;
283  }
284  } while ( gReject < G4UniformRand() );
285 
286  if (G4int(2*G4UniformRand())) epsilon = (1. - epsilon); // Extención de Epsilon hasta 1.
287 
288  G4double logepsilon = log(epsilon);
289  G4double deltaP_R1 = 1. + (a + b*logepsilon + c/logepsilon + d*pow(logepsilon,2.) + e/pow(logepsilon,2.) +
290  f*pow(logepsilon,3.) + gLocal/pow(logepsilon,3.) + h*pow(logepsilon,4.) + i/pow(logepsilon,4.) +
291  j*pow(logepsilon,5.) + k/pow(logepsilon,5.))/100.;
292  G4double deltaP_R2 = 1.+((aa + cc*logepsilon + ee*pow(logepsilon,2.) + gg*pow(logepsilon,3.) + ii*pow(logepsilon,4.))
293  / (1. + bb*logepsilon + dd*pow(logepsilon,2.) + ff*pow(logepsilon,3.) + hh*pow(logepsilon,4.)
294  + jj*pow(logepsilon,5.) ))/100.;
295 
296  if (epsilon <= 0.5)
297  {
298  Rechazo = deltaP_R1/NormaRC;
299  }
300  else
301  {
302  Rechazo = deltaP_R2/NormaRC;
303  }
304  G4cout << Rechazo << " " << NormaRC << " " << epsilon << G4endl;
305  } while (Rechazo < G4UniformRand() );
306 
307  electronTotEnergy = (1. - epsilon) * photonEnergy;
308  positronTotEnergy = epsilon * photonEnergy;
309 
310  } // End of epsilon sampling
311 
312  // Fix charges randomly
313 
314  // Scattered electron (positron) angles. ( Z - axis along the parent photon)
315  // Universal distribution suggested by L. Urban (Geant3 manual (1993) Phys211),
316  // derived from Tsai distribution (Rev. Mod. Phys. 49, 421 (1977)
317 
318  G4double u;
319  const G4double a1 = 0.625;
320  G4double a2 = 3. * a1;
321  // G4double d = 27. ;
322 
323  // if (9. / (9. + d) > G4UniformRand())
324  if (0.25 > G4UniformRand())
325  {
326  u = - std::log(G4UniformRand() * G4UniformRand()) / a1 ;
327  }
328  else
329  {
330  u = - std::log(G4UniformRand() * G4UniformRand()) / a2 ;
331  }
332 
333  G4double thetaEle = u*electron_mass_c2/electronTotEnergy;
334  G4double thetaPos = u*electron_mass_c2/positronTotEnergy;
335  G4double phi = twopi * G4UniformRand();
336 
337  G4double dxEle= std::sin(thetaEle)*std::cos(phi),dyEle= std::sin(thetaEle)*std::sin(phi),dzEle=std::cos(thetaEle);
338  G4double dxPos=-std::sin(thetaPos)*std::cos(phi),dyPos=-std::sin(thetaPos)*std::sin(phi),dzPos=std::cos(thetaPos);
339 
340 
341  // Kinematics of the created pair:
342  // the electron and positron are assumed to have a symetric angular
343  // distribution with respect to the Z axis along the parent photon
344 
345  G4double electronKineEnergy = std::max(0.,electronTotEnergy - electron_mass_c2) ;
346 
347  // SI - The range test has been removed wrt original G4LowEnergyGammaconversion class
348 
349  G4ThreeVector electronDirection (dxEle, dyEle, dzEle);
350  electronDirection.rotateUz(photonDirection);
351 
353  electronDirection,
354  electronKineEnergy);
355 
356  // The e+ is always created (even with kinetic energy = 0) for further annihilation
357  G4double positronKineEnergy = std::max(0.,positronTotEnergy - electron_mass_c2) ;
358 
359  // SI - The range test has been removed wrt original G4LowEnergyGammaconversion class
360 
361  G4ThreeVector positronDirection (dxPos, dyPos, dzPos);
362  positronDirection.rotateUz(photonDirection);
363 
364  // Create G4DynamicParticle object for the particle2
366  positronDirection,
367  positronKineEnergy);
368  // Fill output vector
369 // G4cout << "Cree el e+ " << epsilon << G4endl;
370  fvect->push_back(particle1);
371  fvect->push_back(particle2);
372 
373  if (HardPhotonEnergy > 0.)
374  {
375  G4double thetaHardPhoton = u*electron_mass_c2/HardPhotonEnergy;
376  phi = twopi * G4UniformRand();
377  G4double dxHardP= std::sin(thetaHardPhoton)*std::cos(phi);
378  G4double dyHardP= std::sin(thetaHardPhoton)*std::sin(phi);
379  G4double dzHardP =std::cos(thetaHardPhoton);
380 
381  G4ThreeVector hardPhotonDirection (dxHardP, dyHardP, dzHardP);
382  hardPhotonDirection.rotateUz(photonDirection);
384  hardPhotonDirection,
385  HardPhotonEnergy);
386  fvect->push_back(particle3);
387  }
388 
389 
390  // kill incident photon
393 
394 }
G4double GetKineticEnergy() const
G4double GetfCoulomb() const
Definition: G4Element.hh:190
G4ParticleDefinition * GetDefinition() const
int G4int
Definition: G4Types.hh:78
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
float electron_mass_c2
Definition: hepunit.py:274
static G4Gamma * Gamma()
Definition: G4Gamma.cc:86
G4double GetlogZ3() const
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:198
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
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
G4double GetZ3() const
const G4Element * SelectRandomAtom(const G4MaterialCutsCouple *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.hh:510

Field Documentation

G4ParticleChangeForGamma* G4LivermoreGammaConversionModelRC::fParticleChange
protected

Definition at line 72 of file G4LivermoreGammaConversionModelRC.hh.

Referenced by Initialise(), and SampleSecondaries().


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