G4MollerBhabhaModel Class Reference

#include <G4MollerBhabhaModel.hh>

Inheritance diagram for G4MollerBhabhaModel:

G4VEmModel G4PolarizedMollerBhabhaModel

Public Member Functions

 G4MollerBhabhaModel (const G4ParticleDefinition *p=0, const G4String &nam="MollerBhabha")
virtual ~G4MollerBhabhaModel ()
virtual void Initialise (const G4ParticleDefinition *, const G4DataVector &)
virtual G4double ComputeCrossSectionPerElectron (const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
virtual G4double ComputeCrossSectionPerAtom (const G4ParticleDefinition *, G4double kineticEnergy, G4double Z, G4double A, G4double cutEnergy, G4double maxEnergy)
virtual G4double CrossSectionPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy, G4double maxEnergy)
virtual G4double ComputeDEDXPerVolume (const G4Material *, const G4ParticleDefinition *, G4double kineticEnergy, G4double cutEnergy)
virtual void SampleSecondaries (std::vector< G4DynamicParticle * > *, const G4MaterialCutsCouple *, const G4DynamicParticle *, G4double tmin, G4double maxEnergy)

Protected Member Functions

virtual G4double MaxSecondaryEnergy (const G4ParticleDefinition *, G4double kinEnergy)
void SetParticle (const G4ParticleDefinition *p)

Protected Attributes

const G4ParticleDefinitionparticle
G4ParticleDefinitiontheElectron
G4ParticleChangeForLossfParticleChange
G4bool isElectron
G4double twoln10
G4double lowLimit

Detailed Description

Definition at line 63 of file G4MollerBhabhaModel.hh.


Constructor & Destructor Documentation

G4MollerBhabhaModel::G4MollerBhabhaModel ( const G4ParticleDefinition p = 0,
const G4String nam = "MollerBhabha" 
)

Definition at line 74 of file G4MollerBhabhaModel.cc.

References G4Electron::Electron(), fParticleChange, SetParticle(), and theElectron.

00076   : G4VEmModel(nam),
00077     particle(0),
00078     isElectron(true),
00079     twoln10(2.0*log(10.0)),
00080     lowLimit(0.02*keV),
00081     isInitialised(false)
00082 {
00083   theElectron = G4Electron::Electron();
00084   if(p) { SetParticle(p); }
00085   fParticleChange = 0;
00086 }

G4MollerBhabhaModel::~G4MollerBhabhaModel (  )  [virtual]

Definition at line 90 of file G4MollerBhabhaModel.cc.

00091 {}


Member Function Documentation

G4double G4MollerBhabhaModel::ComputeCrossSectionPerAtom ( const G4ParticleDefinition ,
G4double  kineticEnergy,
G4double  Z,
G4double  A,
G4double  cutEnergy,
G4double  maxEnergy 
) [virtual]

Reimplemented from G4VEmModel.

Definition at line 172 of file G4MollerBhabhaModel.cc.

References ComputeCrossSectionPerElectron().

00178 {
00179   return Z*ComputeCrossSectionPerElectron(p,kineticEnergy,cutEnergy,maxEnergy);
00180 }

G4double G4MollerBhabhaModel::ComputeCrossSectionPerElectron ( const G4ParticleDefinition ,
G4double  kineticEnergy,
G4double  cutEnergy,
G4double  maxEnergy 
) [virtual]

Reimplemented in G4PolarizedMollerBhabhaModel.

Definition at line 119 of file G4MollerBhabhaModel.cc.

References G4InuclParticleNames::gam, isElectron, MaxSecondaryEnergy(), particle, and SetParticle().

Referenced by ComputeCrossSectionPerAtom(), G4PolarizedMollerBhabhaModel::ComputeCrossSectionPerElectron(), and CrossSectionPerVolume().

00123 {
00124   if(!particle) { SetParticle(p); }
00125 
00126   G4double cross = 0.0;
00127   G4double tmax = MaxSecondaryEnergy(p, kineticEnergy);
00128   tmax = std::min(maxEnergy, tmax);
00129 
00130   if(cutEnergy < tmax) {
00131 
00132     G4double xmin  = cutEnergy/kineticEnergy;
00133     G4double xmax  = tmax/kineticEnergy;
00134     G4double tau   = kineticEnergy/electron_mass_c2;
00135     G4double gam   = tau + 1.0;
00136     G4double gamma2= gam*gam;
00137     G4double beta2 = tau*(tau + 2)/gamma2;
00138 
00139     //Moller (e-e-) scattering
00140     if (isElectron) {
00141 
00142       G4double gg = (2.0*gam - 1.0)/gamma2;
00143       cross = ((xmax - xmin)*(1.0 - gg + 1.0/(xmin*xmax)
00144                               + 1.0/((1.0-xmin)*(1.0 - xmax)))
00145             - gg*log( xmax*(1.0 - xmin)/(xmin*(1.0 - xmax)) ) ) / beta2;
00146 
00147     //Bhabha (e+e-) scattering
00148     } else {
00149 
00150       G4double y   = 1.0/(1.0 + gam);
00151       G4double y2  = y*y;
00152       G4double y12 = 1.0 - 2.0*y;
00153       G4double b1  = 2.0 - y2;
00154       G4double b2  = y12*(3.0 + y2);
00155       G4double y122= y12*y12;
00156       G4double b4  = y122*y12;
00157       G4double b3  = b4 + y122;
00158 
00159       cross = (xmax - xmin)*(1.0/(beta2*xmin*xmax) + b2
00160             - 0.5*b3*(xmin + xmax)
00161             + b4*(xmin*xmin + xmin*xmax + xmax*xmax)/3.0)
00162             - b1*log(xmax/xmin);
00163     }
00164 
00165     cross *= twopi_mc2_rcl2/kineticEnergy;
00166   }
00167   return cross;
00168 }

G4double G4MollerBhabhaModel::ComputeDEDXPerVolume ( const G4Material ,
const G4ParticleDefinition ,
G4double  kineticEnergy,
G4double  cutEnergy 
) [virtual]

Reimplemented from G4VEmModel.

Definition at line 210 of file G4MollerBhabhaModel.cc.

References G4IonisParamMat::DensityCorrection(), G4InuclParticleNames::gam, G4Material::GetElectronDensity(), G4Material::GetIonisation(), G4IonisParamMat::GetMeanExcitationEnergy(), G4Material::GetTotNbOfAtomsPerVolume(), isElectron, MaxSecondaryEnergy(), particle, SetParticle(), and twoln10.

00215 {
00216   if(!particle) { SetParticle(p); }
00217   // calculate the dE/dx due to the ionization by Seltzer-Berger formula
00218   // checl low-energy limit
00219   G4double electronDensity = material->GetElectronDensity();
00220   
00221   G4double Zeff  = electronDensity/material->GetTotNbOfAtomsPerVolume();
00222   G4double th    = 0.25*sqrt(Zeff)*keV;
00223   //  G4double cut;  
00224   // if(isElectron) { cut  = std::max(th*0.5, cutEnergy); }
00225   // else           { cut  = std::max(th, cutEnergy); }
00226   G4double tkin  = kineticEnergy;
00227   if (kineticEnergy < th) { tkin = th; }
00228  
00229   G4double tau   = tkin/electron_mass_c2;
00230   G4double gam   = tau + 1.0;
00231   G4double gamma2= gam*gam;
00232   G4double bg2   = tau*(tau + 2);
00233   G4double beta2 = bg2/gamma2;
00234 
00235   G4double eexc  = material->GetIonisation()->GetMeanExcitationEnergy();
00236   eexc          /= electron_mass_c2;
00237   G4double eexc2 = eexc*eexc; 
00238 
00239   G4double d = std::min(cut, MaxSecondaryEnergy(p, tkin))/electron_mass_c2;
00240   G4double dedx;
00241 
00242   // electron
00243   if (isElectron) {
00244 
00245     dedx = log(2.0*(tau + 2.0)/eexc2) - 1.0 - beta2
00246          + log((tau-d)*d) + tau/(tau-d)
00247          + (0.5*d*d + (2.0*tau + 1.)*log(1. - d/tau))/gamma2;
00248    
00249   //positron
00250   } else {
00251 
00252     G4double d2 = d*d*0.5;
00253     G4double d3 = d2*d/1.5;
00254     G4double d4 = d3*d*0.75;
00255     G4double y  = 1.0/(1.0 + gam);
00256     dedx = log(2.0*(tau + 2.0)/eexc2) + log(tau*d)
00257          - beta2*(tau + 2.0*d - y*(3.0*d2 
00258          + y*(d - d3 + y*(d2 - tau*d3 + d4))))/tau;
00259   } 
00260 
00261   //density correction 
00262   G4double x = log(bg2)/twoln10;
00263   dedx -= material->GetIonisation()->DensityCorrection(x); 
00264 
00265   // now you can compute the total ionization loss
00266   dedx *= twopi_mc2_rcl2*electronDensity/beta2;
00267   if (dedx < 0.0) { dedx = 0.0; }
00268  
00269   // lowenergy extrapolation
00270  
00271   if (kineticEnergy < th) {
00272     x = kineticEnergy/th;
00273     if(x > 0.25) { dedx /= sqrt(x); }
00274     else { dedx *= 1.4*sqrt(x)/(0.1 + x); }
00275   }
00276   return dedx;
00277 }

G4double G4MollerBhabhaModel::CrossSectionPerVolume ( const G4Material ,
const G4ParticleDefinition ,
G4double  kineticEnergy,
G4double  cutEnergy,
G4double  maxEnergy 
) [virtual]

Reimplemented from G4VEmModel.

Definition at line 184 of file G4MollerBhabhaModel.cc.

References ComputeCrossSectionPerElectron(), and G4Material::GetElectronDensity().

00190 {
00191   G4double eDensity = material->GetElectronDensity();
00192   return eDensity*ComputeCrossSectionPerElectron(p,kinEnergy,cutEnergy,maxEnergy);
00193   /*
00194   G4double Zeff     = eDensity/material->GetTotNbOfAtomsPerVolume();
00195   G4double th       = 0.25*sqrt(Zeff)*keV;
00196   G4double cut;  
00197   if(isElectron) { cut  = std::max(th*0.5, cutEnergy); }
00198   else           { cut  = std::max(th, cutEnergy); }
00199   G4double res = 0.0;
00200   // below this threshold no bremsstrahlung
00201   if (kinEnergy > th) { 
00202     res = eDensity*ComputeCrossSectionPerElectron(p,kinEnergy,cut,maxEnergy);
00203   }
00204   return res; 
00205   */
00206 }

void G4MollerBhabhaModel::Initialise ( const G4ParticleDefinition ,
const G4DataVector  
) [virtual]

Implements G4VEmModel.

Definition at line 105 of file G4MollerBhabhaModel.cc.

References fParticleChange, G4VEmModel::GetParticleChangeForLoss(), particle, and SetParticle().

00107 {
00108   if(!particle) { SetParticle(p); }
00109 
00110   if(isInitialised) { return; }
00111 
00112   isInitialised = true;
00113   fParticleChange = GetParticleChangeForLoss();
00114 }

G4double G4MollerBhabhaModel::MaxSecondaryEnergy ( const G4ParticleDefinition ,
G4double  kinEnergy 
) [protected, virtual]

Reimplemented from G4VEmModel.

Definition at line 95 of file G4MollerBhabhaModel.cc.

References isElectron.

Referenced by G4PolarizedMollerBhabhaModel::ComputeCrossSectionPerElectron(), ComputeCrossSectionPerElectron(), and ComputeDEDXPerVolume().

00097 {
00098   G4double tmax = kinEnergy;
00099   if(isElectron) { tmax *= 0.5; }
00100   return tmax;
00101 }

void G4MollerBhabhaModel::SampleSecondaries ( std::vector< G4DynamicParticle * > *  ,
const G4MaterialCutsCouple ,
const G4DynamicParticle ,
G4double  tmin,
G4double  maxEnergy 
) [virtual]

Implements G4VEmModel.

Reimplemented in G4PolarizedMollerBhabhaModel.

Definition at line 281 of file G4MollerBhabhaModel.cc.

References fParticleChange, G4UniformRand, G4InuclParticleNames::gam, G4DynamicParticle::GetKineticEnergy(), G4DynamicParticle::GetMomentumDirection(), isElectron, G4ParticleChangeForLoss::SetProposedKineticEnergy(), G4ParticleChangeForLoss::SetProposedMomentumDirection(), and theElectron.

00286 {
00287   G4double kineticEnergy = dp->GetKineticEnergy();
00288   //const G4Material* mat = couple->GetMaterial();
00289   //G4double Zeff = mat->GetElectronDensity()/mat->GetTotNbOfAtomsPerVolume();
00290   // G4double th   = 0.25*sqrt(Zeff)*keV;
00291   G4double tmax;
00292   G4double tmin = cutEnergy;  
00293   if(isElectron) { 
00294     tmax = 0.5*kineticEnergy; 
00295   } else {
00296     tmax = kineticEnergy; 
00297   }
00298   if(maxEnergy < tmax) { tmax = maxEnergy; }
00299   if(tmin >= tmax) { return; }
00300 
00301   G4double energy = kineticEnergy + electron_mass_c2;
00302   G4double totalMomentum = sqrt(kineticEnergy*(energy + electron_mass_c2));
00303   G4double xmin   = tmin/kineticEnergy;
00304   G4double xmax   = tmax/kineticEnergy;
00305   G4double gam    = energy/electron_mass_c2;
00306   G4double gamma2 = gam*gam;
00307   G4double beta2  = 1.0 - 1.0/gamma2;
00308   G4double x, z, q, grej;
00309 
00310   G4ThreeVector direction = dp->GetMomentumDirection();
00311 
00312   //Moller (e-e-) scattering
00313   if (isElectron) {
00314 
00315     G4double gg = (2.0*gam - 1.0)/gamma2;
00316     G4double y = 1.0 - xmax;
00317     grej = 1.0 - gg*xmax + xmax*xmax*(1.0 - gg + (1.0 - gg*y)/(y*y));
00318 
00319     do {
00320       q = G4UniformRand();
00321       x = xmin*xmax/(xmin*(1.0 - q) + xmax*q);
00322       y = 1.0 - x;
00323       z = 1.0 - gg*x + x*x*(1.0 - gg + (1.0 - gg*y)/(y*y));
00324       /*
00325       if(z > grej) {
00326         G4cout << "G4MollerBhabhaModel::SampleSecondary Warning! "
00327                << "Majorant " << grej << " < "
00328                << z << " for x= " << x
00329                << " e-e- scattering"
00330                << G4endl;
00331       }
00332       */
00333     } while(grej * G4UniformRand() > z);
00334 
00335   //Bhabha (e+e-) scattering
00336   } else {
00337 
00338     G4double y   = 1.0/(1.0 + gam);
00339     G4double y2  = y*y;
00340     G4double y12 = 1.0 - 2.0*y;
00341     G4double b1  = 2.0 - y2;
00342     G4double b2  = y12*(3.0 + y2);
00343     G4double y122= y12*y12;
00344     G4double b4  = y122*y12;
00345     G4double b3  = b4 + y122;
00346 
00347     y    = xmax*xmax;
00348     grej = 1.0 + (y*y*b4 - xmin*xmin*xmin*b3 + y*b2 - xmin*b1)*beta2; 
00349     do {
00350       q = G4UniformRand();
00351       x = xmin*xmax/(xmin*(1.0 - q) + xmax*q);
00352       y = x*x;
00353       z = 1.0 + (y*y*b4 - x*y*b3 + y*b2 - x*b1)*beta2; 
00354       /*
00355       if(z > grej) {
00356         G4cout << "G4MollerBhabhaModel::SampleSecondary Warning! "
00357                << "Majorant " << grej << " < "
00358                << z << " for x= " << x
00359                << " e+e- scattering"
00360                << G4endl;
00361       }
00362       */
00363     } while(grej * G4UniformRand() > z);
00364   }
00365 
00366   G4double deltaKinEnergy = x * kineticEnergy;
00367 
00368   G4double deltaMomentum =
00369            sqrt(deltaKinEnergy * (deltaKinEnergy + 2.0*electron_mass_c2));
00370   G4double cost = deltaKinEnergy * (energy + electron_mass_c2) /
00371                                    (deltaMomentum * totalMomentum);
00372   G4double sint = (1.0 - cost)*(1. + cost);
00373   if(sint > 0.0) { sint = sqrt(sint); }
00374   else { sint = 0.0; }
00375 
00376   G4double phi = twopi * G4UniformRand() ;
00377 
00378   G4ThreeVector deltaDirection(sint*cos(phi),sint*sin(phi), cost) ;
00379   deltaDirection.rotateUz(direction);
00380 
00381   // primary change
00382   kineticEnergy -= deltaKinEnergy;
00383   fParticleChange->SetProposedKineticEnergy(kineticEnergy);
00384 
00385   G4ThreeVector dir = totalMomentum*direction - deltaMomentum*deltaDirection;
00386   direction = dir.unit();
00387   fParticleChange->SetProposedMomentumDirection(direction);
00388 
00389   // create G4DynamicParticle object for delta ray
00390   G4DynamicParticle* delta = new G4DynamicParticle(theElectron,
00391                                                    deltaDirection,deltaKinEnergy);
00392   vdp->push_back(delta);
00393 }

void G4MollerBhabhaModel::SetParticle ( const G4ParticleDefinition p  )  [inline, protected]

Definition at line 132 of file G4MollerBhabhaModel.hh.

References isElectron, particle, and theElectron.

Referenced by ComputeCrossSectionPerElectron(), ComputeDEDXPerVolume(), G4MollerBhabhaModel(), and Initialise().

00133 {
00134   particle = p;
00135   if(p != theElectron) { isElectron = false; }
00136 }


Field Documentation

G4ParticleChangeForLoss* G4MollerBhabhaModel::fParticleChange [protected]

Definition at line 114 of file G4MollerBhabhaModel.hh.

Referenced by G4MollerBhabhaModel(), Initialise(), G4PolarizedMollerBhabhaModel::SampleSecondaries(), and SampleSecondaries().

G4bool G4MollerBhabhaModel::isElectron [protected]

Definition at line 116 of file G4MollerBhabhaModel.hh.

Referenced by ComputeCrossSectionPerElectron(), ComputeDEDXPerVolume(), G4PolarizedMollerBhabhaModel::G4PolarizedMollerBhabhaModel(), MaxSecondaryEnergy(), G4PolarizedMollerBhabhaModel::SampleSecondaries(), SampleSecondaries(), and SetParticle().

G4double G4MollerBhabhaModel::lowLimit [protected]

Reimplemented from G4VEmModel.

Definition at line 118 of file G4MollerBhabhaModel.hh.

const G4ParticleDefinition* G4MollerBhabhaModel::particle [protected]

Definition at line 112 of file G4MollerBhabhaModel.hh.

Referenced by ComputeCrossSectionPerElectron(), ComputeDEDXPerVolume(), Initialise(), and SetParticle().

G4ParticleDefinition* G4MollerBhabhaModel::theElectron [protected]

Definition at line 113 of file G4MollerBhabhaModel.hh.

Referenced by G4MollerBhabhaModel(), G4PolarizedMollerBhabhaModel::G4PolarizedMollerBhabhaModel(), G4PolarizedMollerBhabhaModel::SampleSecondaries(), SampleSecondaries(), and SetParticle().

G4double G4MollerBhabhaModel::twoln10 [protected]

Definition at line 117 of file G4MollerBhabhaModel.hh.

Referenced by ComputeDEDXPerVolume().


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:52:30 2013 for Geant4 by  doxygen 1.4.7