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

#include <G4LivermorePolarizedGammaConversionModel.hh>

Inheritance diagram for G4LivermorePolarizedGammaConversionModel:
G4VEmModel

Public Member Functions

 G4LivermorePolarizedGammaConversionModel (const G4ParticleDefinition *p=0, const G4String &nam="LivermorePolarizedGammaConversion")
 
virtual ~G4LivermorePolarizedGammaConversionModel ()
 
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 45 of file G4LivermorePolarizedGammaConversionModel.hh.

Constructor & Destructor Documentation

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

Definition at line 42 of file G4LivermorePolarizedGammaConversionModel.cc.

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

44  :G4VEmModel(nam),fParticleChange(0),
45  isInitialised(false),meanFreePathTable(0),crossSectionHandler(0)
46 {
47  lowEnergyLimit = 2*electron_mass_c2;
48  highEnergyLimit = 100 * GeV;
49  SetLowEnergyLimit(lowEnergyLimit);
50  SetHighEnergyLimit(highEnergyLimit);
51  smallEnergy = 2.*MeV;
52 
53  Phi=0.;
54  Psi=0.;
55 
56  verboseLevel= 0;
57  // Verbosity scale:
58  // 0 = nothing
59  // 1 = warning for energy non-conservation
60  // 2 = details of energy budget
61  // 3 = calculation of cross sections, file openings, samping of atoms
62  // 4 = entering in methods
63 
64  if(verboseLevel > 0) {
65  G4cout << "Livermore Polarized GammaConversion is constructed " << G4endl
66  << "Energy range: "
67  << lowEnergyLimit / keV << " keV - "
68  << highEnergyLimit / GeV << " GeV"
69  << G4endl;
70  }
71 
72  crossSectionHandler = new G4CrossSectionHandler();
73  crossSectionHandler->Initialise(0,lowEnergyLimit,highEnergyLimit,400);
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
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)
#define G4endl
Definition: G4ios.hh:61
void SetLowEnergyLimit(G4double)
Definition: G4VEmModel.hh:690
G4LivermorePolarizedGammaConversionModel::~G4LivermorePolarizedGammaConversionModel ( )
virtual

Definition at line 78 of file G4LivermorePolarizedGammaConversionModel.cc.

79 {
80  delete crossSectionHandler;
81 }

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

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

153 {
154  if (verboseLevel > 3) {
155  G4cout << "G4LivermorePolarizedGammaConversionModel::ComputeCrossSectionPerAtom()"
156  << G4endl;
157  }
158  if(Z < 0.9 || GammaEnergy <= lowEnergyLimit) { return 0.0; }
159  G4double cs = crossSectionHandler->FindValue(G4int(Z), GammaEnergy);
160  return cs;
161 }
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 G4LivermorePolarizedGammaConversionModel::GetMeanFreePath ( const G4Track aTrack,
G4double  previousStepSize,
G4ForceCondition condition 
)
protected
void G4LivermorePolarizedGammaConversionModel::Initialise ( const G4ParticleDefinition particle,
const G4DataVector cuts 
)
virtual

Implements G4VEmModel.

Definition at line 87 of file G4LivermorePolarizedGammaConversionModel.cc.

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

89 {
90  if (verboseLevel > 3)
91  G4cout << "Calling G4LivermorePolarizedGammaConversionModel::Initialise()" << G4endl;
92 
93  if (crossSectionHandler)
94  {
95  crossSectionHandler->Clear();
96  delete crossSectionHandler;
97  }
98 
99  // Energy limits
100  /*
101  // V.Ivanchenko: this was meanless check
102  if (LowEnergyLimit() < lowEnergyLimit)
103  {
104  G4cout << "G4LivermorePolarizedGammaConversionModel: low energy limit increased from " <<
105  LowEnergyLimit()/eV << " eV to " << lowEnergyLimit << " eV" << G4endl;
106  // SetLowEnergyLimit(lowEnergyLimit);
107  }
108  */
109  if (HighEnergyLimit() > highEnergyLimit)
110  {
111  G4cout << "G4LivermorePolarizedGammaConversionModel: high energy limit decreased from " <<
112  HighEnergyLimit()/GeV << " GeV to " << highEnergyLimit << " GeV" << G4endl;
113  // V.Ivanchenko: this is forbidden
114  // SetHighEnergyLimit(highEnergyLimit);
115  }
116 
117  // Reading of data files - all materials are read
118 
119  crossSectionHandler = new G4CrossSectionHandler;
120  crossSectionHandler->Clear();
121  G4String crossSectionFile = "pair/pp-cs-";
122  crossSectionHandler->LoadData(crossSectionFile);
123 
124  //
125  if (verboseLevel > 2) {
126  G4cout << "Loaded cross section files for Livermore Polarized GammaConversion model"
127  << G4endl;
128  }
129  InitialiseElementSelectors(particle,cuts);
130 
131  if(verboseLevel > 0) {
132  G4cout << "Livermore Polarized GammaConversion model is initialized " << G4endl
133  << "Energy range: "
134  << LowEnergyLimit() / keV << " keV - "
135  << HighEnergyLimit() / GeV << " GeV"
136  << G4endl;
137  }
138 
139  //
140  if(!isInitialised) {
141  isInitialised = true;
143  }
144 }
G4double LowEnergyLimit() const
Definition: G4VEmModel.hh:599
void InitialiseElementSelectors(const G4ParticleDefinition *, const G4DataVector &)
Definition: G4VEmModel.cc:135
G4double HighEnergyLimit() const
Definition: G4VEmModel.hh:592
G4GLOB_DLL std::ostream G4cout
void LoadData(const G4String &dataFile)
#define G4endl
Definition: G4ios.hh:61
G4ParticleChangeForGamma * GetParticleChangeForGamma()
Definition: G4VEmModel.cc:121
void G4LivermorePolarizedGammaConversionModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple couple,
const G4DynamicParticle aDynamicGamma,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Implements G4VEmModel.

Definition at line 166 of file G4LivermorePolarizedGammaConversionModel.cc.

References choice(), G4Electron::Electron(), python.hepunit::electron_mass_c2, fParticleChange, fStopAndKill, G4cout, G4endl, G4UniformRand, G4DynamicParticle::GetDefinition(), G4Element::GetfCoulomb(), G4Element::GetIonisation(), G4DynamicParticle::GetKineticEnergy(), G4IonisParamElm::GetlogZ3(), G4DynamicParticle::GetMomentumDirection(), G4DynamicParticle::GetPolarization(), G4IonisParamElm::GetZ3(), CLHEP::Hep3Vector::howOrthogonal(), CLHEP::Hep3Vector::isOrthogonal(), CLHEP::Hep3Vector::mag(), G4INCL::Math::max(), python.hepunit::MeV, G4INCL::Math::min(), G4Positron::Positron(), G4ParticleChangeForGamma::ProposeMomentumDirection(), G4VParticleChange::ProposeTrackStatus(), G4VEmModel::SelectRandomAtom(), and G4ParticleChangeForGamma::SetProposedKineticEnergy().

171 {
172 
173  // Fluorescence generated according to:
174  // J. Stepanek ,"A program to determine the radiation spectra due to a single atomic
175  // subshell ionisation by a particle or due to deexcitation or decay of radionuclides",
176  // Comp. Phys. Comm. 1206 pp 1-1-9 (1997)
177 
178  if (verboseLevel > 3)
179  G4cout << "Calling SampleSecondaries() of G4LivermorePolarizedGammaConversionModel" << G4endl;
180 
181  G4double photonEnergy = aDynamicGamma->GetKineticEnergy();
182  // Within energy limit?
183 
184  if(photonEnergy <= lowEnergyLimit)
185  {
188  return;
189  }
190 
191 
192  G4ThreeVector gammaPolarization0 = aDynamicGamma->GetPolarization();
193  G4ThreeVector gammaDirection0 = aDynamicGamma->GetMomentumDirection();
194 
195  // Make sure that the polarization vector is perpendicular to the
196  // gamma direction. If not
197 
198  if(!(gammaPolarization0.isOrthogonal(gammaDirection0, 1e-6))||(gammaPolarization0.mag()==0))
199  { // only for testing now
200  gammaPolarization0 = GetRandomPolarization(gammaDirection0);
201  }
202  else
203  {
204  if ( gammaPolarization0.howOrthogonal(gammaDirection0) != 0)
205  {
206  gammaPolarization0 = GetPerpendicularPolarization(gammaDirection0, gammaPolarization0);
207  }
208  }
209 
210  // End of Protection
211 
212 
213  G4double epsilon ;
214  G4double epsilon0Local = electron_mass_c2 / photonEnergy ;
215 
216  // Do it fast if photon energy < 2. MeV
217 
218  if (photonEnergy < smallEnergy )
219  {
220  epsilon = epsilon0Local + (0.5 - epsilon0Local) * G4UniformRand();
221  }
222  else
223  {
224 
225  // Select randomly one element in the current material
226 
227  // G4int Z = crossSectionHandler->SelectRandomAtom(couple,photonEnergy);
228  //const G4Element* element = crossSectionHandler->SelectRandomElement(couple,photonEnergy);
229 
230  const G4ParticleDefinition* particle = aDynamicGamma->GetDefinition();
231  const G4Element* element = SelectRandomAtom(couple,particle,photonEnergy);
232 
233  /*
234  if (element == 0)
235  {
236  G4cout << "G4LivermorePolarizedGammaConversionModel::PostStepDoIt - element = 0" << G4endl;
237  }
238  */
239 
240  G4IonisParamElm* ionisation = element->GetIonisation();
241 
242  /*
243  if (ionisation == 0)
244  {
245  G4cout << "G4LivermorePolarizedGammaConversionModel::PostStepDoIt - ionisation = 0" << G4endl;
246  }
247  */
248 
249 
250  // Extract Coulomb factor for this Element
251 
252  G4double fZ = 8. * (ionisation->GetlogZ3());
253  if (photonEnergy > 50. * MeV) fZ += 8. * (element->GetfCoulomb());
254 
255  // Limits of the screening variable
256  G4double screenFactor = 136. * epsilon0Local / (element->GetIonisation()->GetZ3()) ;
257  G4double screenMax = exp ((42.24 - fZ)/8.368) - 0.952 ;
258  G4double screenMin = std::min(4.*screenFactor,screenMax) ;
259 
260  // Limits of the energy sampling
261  G4double epsilon1 = 0.5 - 0.5 * sqrt(1. - screenMin / screenMax) ;
262  G4double epsilonMin = std::max(epsilon0Local,epsilon1);
263  G4double epsilonRange = 0.5 - epsilonMin ;
264 
265  // Sample the energy rate of the created electron (or positron)
266  G4double screen;
267  G4double gReject ;
268 
269  G4double f10 = ScreenFunction1(screenMin) - fZ;
270  G4double f20 = ScreenFunction2(screenMin) - fZ;
271  G4double normF1 = std::max(f10 * epsilonRange * epsilonRange,0.);
272  G4double normF2 = std::max(1.5 * f20,0.);
273 
274  do {
275  if (normF1 / (normF1 + normF2) > G4UniformRand() )
276  {
277  epsilon = 0.5 - epsilonRange * pow(G4UniformRand(), 0.3333) ;
278  screen = screenFactor / (epsilon * (1. - epsilon));
279  gReject = (ScreenFunction1(screen) - fZ) / f10 ;
280  }
281  else
282  {
283  epsilon = epsilonMin + epsilonRange * G4UniformRand();
284  screen = screenFactor / (epsilon * (1 - epsilon));
285  gReject = (ScreenFunction2(screen) - fZ) / f20 ;
286 
287 
288  }
289  } while ( gReject < G4UniformRand() );
290 
291  } // End of epsilon sampling
292 
293  // Fix charges randomly
294 
295  G4double electronTotEnergy;
296  G4double positronTotEnergy;
297 
298 
299  if (G4int(2*G4UniformRand()))
300  {
301  electronTotEnergy = (1. - epsilon) * photonEnergy;
302  positronTotEnergy = epsilon * photonEnergy;
303  }
304  else
305  {
306  positronTotEnergy = (1. - epsilon) * photonEnergy;
307  electronTotEnergy = epsilon * photonEnergy;
308  }
309 
310  // Scattered electron (positron) angles. ( Z - axis along the parent photon)
311  // Universal distribution suggested by L. Urban (Geant3 manual (1993) Phys211),
312  // derived from Tsai distribution (Rev. Mod. Phys. 49, 421 (1977)
313 
314 /*
315  G4double u;
316  const G4double a1 = 0.625;
317  G4double a2 = 3. * a1;
318 
319  if (0.25 > G4UniformRand())
320  {
321  u = - log(G4UniformRand() * G4UniformRand()) / a1 ;
322  }
323  else
324  {
325  u = - log(G4UniformRand() * G4UniformRand()) / a2 ;
326  }
327 */
328 
329  G4double Ene = electronTotEnergy/electron_mass_c2; // Normalized energy
330 
331  G4double cosTheta = 0.;
332  G4double sinTheta = 0.;
333 
334  SetTheta(&cosTheta,&sinTheta,Ene);
335 
336  // G4double theta = u * electron_mass_c2 / photonEnergy ;
337  // G4double phi = twopi * G4UniformRand() ;
338 
339  G4double phi,psi=0.;
340 
341  //corrected e+ e- angular angular distribution //preliminary!
342 
343  // if(photonEnergy>50*MeV)
344  // {
345  phi = SetPhi(photonEnergy);
346  psi = SetPsi(photonEnergy,phi);
347  // }
348  //else
349  // {
350  //psi = G4UniformRand()*2.*pi;
351  //phi = pi; // coplanar
352  // }
353 
354  Psi = psi;
355  Phi = phi;
356  //G4cout << "PHI " << phi << G4endl;
357  //G4cout << "PSI " << psi << G4endl;
358 
359  G4double phie, phip;
361  choice = G4UniformRand();
362  if (choice <= 0.5)
363  {
364  phie = psi; //azimuthal angle for the electron
365  phip = phie+phi; //azimuthal angle for the positron
366  }
367  else
368  {
369  // opzione 1 phie / phip equivalenti
370 
371  phip = psi; //azimuthal angle for the positron
372  phie = phip + phi; //azimuthal angle for the electron
373  }
374 
375 
376  // Electron Kinematics
377 
378  G4double dirX = sinTheta*cos(phie);
379  G4double dirY = sinTheta*sin(phie);
380  G4double dirZ = cosTheta;
381  G4ThreeVector electronDirection(dirX,dirY,dirZ);
382 
383  // Kinematics of the created pair:
384  // the electron and positron are assumed to have a symetric angular
385  // distribution with respect to the Z axis along the parent photon
386 
387  //G4double localEnergyDeposit = 0. ;
388 
389  G4double electronKineEnergy = std::max(0.,electronTotEnergy - electron_mass_c2) ;
390 
391  SystemOfRefChange(gammaDirection0,electronDirection,
392  gammaPolarization0);
393 
395  electronDirection,
396  electronKineEnergy);
397 
398  // The e+ is always created (even with kinetic energy = 0) for further annihilation
399 
400  Ene = positronTotEnergy/electron_mass_c2; // Normalized energy
401 
402  cosTheta = 0.;
403  sinTheta = 0.;
404 
405  SetTheta(&cosTheta,&sinTheta,Ene);
406 
407  // Positron Kinematics
408 
409  dirX = sinTheta*cos(phip);
410  dirY = sinTheta*sin(phip);
411  dirZ = cosTheta;
412  G4ThreeVector positronDirection(dirX,dirY,dirZ);
413 
414  G4double positronKineEnergy = std::max(0.,positronTotEnergy - electron_mass_c2) ;
415  SystemOfRefChange(gammaDirection0,positronDirection,
416  gammaPolarization0);
417 
418  // Create G4DynamicParticle object for the particle2
420  positronDirection, positronKineEnergy);
421 
422 
423  fvect->push_back(particle1);
424  fvect->push_back(particle2);
425 
426 
427 
428  // Kill the incident photon
429 
430 
431 
432  // Create lists of pointers to DynamicParticles (photons and electrons)
433  // (Is the electron vector necessary? To be checked)
434  // std::vector<G4DynamicParticle*>* photonVector = 0;
435  //std::vector<G4DynamicParticle*> electronVector;
436 
440 
441 }
G4double GetKineticEnergy() const
G4double GetfCoulomb() const
Definition: G4Element.hh:190
G4ParticleDefinition * GetDefinition() const
bool isOrthogonal(const Hep3Vector &v, double epsilon=tolerance) const
Definition: SpaceVector.cc:237
int G4int
Definition: G4Types.hh:78
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
double howOrthogonal(const Hep3Vector &v) const
Definition: SpaceVector.cc:219
#define G4UniformRand()
Definition: Randomize.hh:87
G4GLOB_DLL std::ostream G4cout
const G4ThreeVector & GetMomentumDirection() const
float electron_mass_c2
Definition: hepunit.py:274
G4double GetlogZ3() const
subroutine choice(MNUM, RR, ICHAN, PROB1, PROB2, PROB3, AMRX, GAMRX, AMRA, GAMRA, AMRB, GAMRB)
Definition: leptonew.f:1817
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
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
double G4double
Definition: G4Types.hh:76
void ProposeTrackStatus(G4TrackStatus status)
double mag() const
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* G4LivermorePolarizedGammaConversionModel::fParticleChange
protected

Definition at line 73 of file G4LivermorePolarizedGammaConversionModel.hh.

Referenced by Initialise(), and SampleSecondaries().


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