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

#include <MyKleinNishinaCompton.hh>

Inheritance diagram for MyKleinNishinaCompton:
G4KleinNishinaCompton G4VEmModel

Public Member Functions

 MyKleinNishinaCompton (DetectorConstruction *, const G4ParticleDefinition *p=0, const G4String &nam="myKlein-Nishina")
 
 ~MyKleinNishinaCompton ()
 
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kinEnergy, G4double cut, G4double emax)
 
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)
 
void SetCSFactor (G4double factor)
 
- Public Member Functions inherited from G4KleinNishinaCompton
 G4KleinNishinaCompton (const G4ParticleDefinition *p=0, const G4String &nam="Klein-Nishina")
 
virtual ~G4KleinNishinaCompton ()
 
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
 
virtual void InitialiseLocal (const G4ParticleDefinition *, G4VEmModel *masterModel)
 
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kinEnergy, G4double Z, G4double A, G4double cut, G4double emax)
 
- Public Member Functions inherited from G4VEmModel
 G4VEmModel (const G4String &nam)
 
virtual ~G4VEmModel ()
 
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 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 Attributes

DetectorConstructionfDetector
 
MyKleinNishinaMessengerfMessenger
 
G4double fCrossSectionFactor
 
- Protected Attributes inherited from G4KleinNishinaCompton
G4ParticleDefinitiontheGamma
 
G4ParticleDefinitiontheElectron
 
G4ParticleChangeForGammafParticleChange
 
G4double lowestGammaEnergy
 
- Protected Attributes inherited from G4VEmModel
G4ElementDatafElementData
 
G4VParticleChangepParticleChange
 
G4PhysicsTablexSectionTable
 
const std::vector< G4double > * theDensityFactor
 
const std::vector< G4int > * theDensityIdx
 
size_t idxTable
 

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 *)
 

Detailed Description

Definition at line 43 of file MyKleinNishinaCompton.hh.

Constructor & Destructor Documentation

MyKleinNishinaCompton::MyKleinNishinaCompton ( DetectorConstruction det,
const G4ParticleDefinition p = 0,
const G4String nam = "myKlein-Nishina" 
)

Definition at line 49 of file MyKleinNishinaCompton.cc.

References fCrossSectionFactor, and fMessenger.

53 {
56 }
DetectorConstruction * fDetector
G4KleinNishinaCompton(const G4ParticleDefinition *p=0, const G4String &nam="Klein-Nishina")
MyKleinNishinaMessenger * fMessenger
MyKleinNishinaCompton::~MyKleinNishinaCompton ( )

Definition at line 60 of file MyKleinNishinaCompton.cc.

References fMessenger.

61 {
62  delete fMessenger;
63 }
MyKleinNishinaMessenger * fMessenger

Member Function Documentation

G4double MyKleinNishinaCompton::CrossSectionPerVolume ( const G4Material mat,
const G4ParticleDefinition part,
G4double  kinEnergy,
G4double  cut,
G4double  emax 
)
virtual

Reimplemented from G4VEmModel.

Definition at line 67 of file MyKleinNishinaCompton.cc.

References G4VEmModel::CrossSectionPerVolume(), and fCrossSectionFactor.

72 {
73  G4double xsection = G4VEmModel::CrossSectionPerVolume(mat,part,GammaEnergy);
74 
75  return xsection*fCrossSectionFactor;
76 }
virtual G4double CrossSectionPerVolume(const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy=0.0, G4double maxEnergy=DBL_MAX)
Definition: G4VEmModel.cc:245
double G4double
Definition: G4Types.hh:76
void MyKleinNishinaCompton::SampleSecondaries ( std::vector< G4DynamicParticle * > *  fvect,
const G4MaterialCutsCouple ,
const G4DynamicParticle aDynamicGamma,
G4double  tmin,
G4double  maxEnergy 
)
virtual

Reimplemented from G4KleinNishinaCompton.

Definition at line 79 of file MyKleinNishinaCompton.cc.

References DBL_MIN, python.hepunit::electron_mass_c2, G4KleinNishinaCompton::fParticleChange, G4UniformRand, G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), G4ParticleChangeForGamma::ProposeMomentumDirection(), CLHEP::Hep3Vector::rotateUz(), G4ParticleChangeForGamma::SetProposedKineticEnergy(), G4KleinNishinaCompton::theElectron, python.hepunit::twopi, and CLHEP::Hep3Vector::unit().

85 {
86  // The scattered gamma energy is sampled according to Klein - Nishina formula.
87  // The random number techniques of Butcher & Messel are used
88  // (Nuc Phys 20(1960),15).
89  // Note : Effects due to binding of atomic electrons are negliged.
90 
91  G4double gamEnergy0 = aDynamicGamma->GetKineticEnergy();
92  G4double E0_m = gamEnergy0 / electron_mass_c2 ;
93 
94  G4ThreeVector gamDirection0 = aDynamicGamma->GetMomentumDirection();
95 
96  //
97  // sample the energy rate of the scattered gamma
98  //
99 
100  G4double epsilon, epsilonsq, onecost, sint2, greject ;
101 
102  G4double eps0 = 1./(1. + 2.*E0_m);
103  G4double eps0sq = eps0*eps0;
104  G4double alpha1 = - log(eps0);
105  G4double alpha2 = 0.5*(1.- eps0sq);
106 
107  do {
108  if ( alpha1/(alpha1+alpha2) > G4UniformRand() ) {
109  epsilon = exp(-alpha1*G4UniformRand()); // eps0**r
110  epsilonsq = epsilon*epsilon;
111 
112  } else {
113  epsilonsq = eps0sq + (1.- eps0sq)*G4UniformRand();
114  epsilon = sqrt(epsilonsq);
115  };
116 
117  onecost = (1.- epsilon)/(epsilon*E0_m);
118  sint2 = onecost*(2.-onecost);
119  greject = 1. - epsilon*sint2/(1.+ epsilonsq);
120 
121  } while (greject < G4UniformRand());
122 
123  //
124  // scattered gamma angles. ( Z - axis along the parent gamma)
125  //
126 
127  G4double cosTeta = 1. - onecost;
128  G4double sinTeta = sqrt (sint2);
129  G4double Phi = twopi * G4UniformRand();
130  G4double dirx = sinTeta*cos(Phi), diry = sinTeta*sin(Phi), dirz = cosTeta;
131 
132  //
133  // update G4VParticleChange for the scattered gamma
134  //
135  // beam regeneration trick : restore incident beam
136 
137  G4ThreeVector gamDirection1 ( dirx,diry,dirz );
138  gamDirection1.rotateUz(gamDirection0);
139  G4double gamEnergy1 = epsilon*gamEnergy0;
142 
143  //
144  // kinematic of the scattered electron
145  //
146 
147  G4double eKinEnergy = gamEnergy0 - gamEnergy1;
148 
149  if(eKinEnergy > DBL_MIN) {
150  G4ThreeVector eDirection
151  = gamEnergy0*gamDirection0 - gamEnergy1*gamDirection1;
152  eDirection = eDirection.unit();
153 
154  // create G4DynamicParticle object for the electron.
155  G4DynamicParticle* dp
156  = new G4DynamicParticle(theElectron,eDirection,eKinEnergy);
157  fvect->push_back(dp);
158  }
159 }
G4ParticleChangeForGamma * fParticleChange
G4double GetKineticEnergy() const
void ProposeMomentumDirection(G4double Px, G4double Py, G4double Pz)
#define G4UniformRand()
Definition: Randomize.hh:87
const G4ThreeVector & GetMomentumDirection() const
float electron_mass_c2
Definition: hepunit.py:274
Hep3Vector unit() const
#define DBL_MIN
Definition: templates.hh:75
void SetProposedKineticEnergy(G4double proposedKinEnergy)
G4ParticleDefinition * theElectron
double G4double
Definition: G4Types.hh:76
void MyKleinNishinaCompton::SetCSFactor ( G4double  factor)
inline

Field Documentation

G4double MyKleinNishinaCompton::fCrossSectionFactor
protected
DetectorConstruction* MyKleinNishinaCompton::fDetector
protected

Definition at line 67 of file MyKleinNishinaCompton.hh.

MyKleinNishinaMessenger* MyKleinNishinaCompton::fMessenger
protected

Definition at line 72 of file MyKleinNishinaCompton.hh.

Referenced by MyKleinNishinaCompton(), and ~MyKleinNishinaCompton().


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