G4VEnergyLoss Class Reference

#include <G4VEnergyLoss.hh>

Inheritance diagram for G4VEnergyLoss:

G4VContinuousDiscreteProcess G4VProcess

Public Member Functions

 G4VEnergyLoss (const G4String &, G4ProcessType aType=fNotDefined)
virtual ~G4VEnergyLoss ()
virtual G4double GetContinuousStepLimit (const G4Track &track, G4double previousStepSize, G4double currentMinimumStep, G4double &currentSafety)=0
virtual G4VParticleChangeAlongStepDoIt (const G4Track &track, const G4Step &Step)=0
virtual G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, G4ForceCondition *condition)=0
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &Step)=0

Static Public Member Functions

static void SetRndmStep (G4bool value)
static void SetEnlossFluc (G4bool value)
static void SetSubSec (G4bool value)
static void SetMinDeltaCutInRange (G4double value)
static void SetStepFunction (G4double c1, G4double c2)

Protected Member Functions

G4double GetLossWithFluct (const G4DynamicParticle *aParticle, const G4MaterialCutsCouple *couple, G4double ChargeSquare, G4double MeanLoss, G4double step)
G4PhysicsTableBuildRangeTable (G4PhysicsTable *theDEDXTable, G4PhysicsTable *theRangeTable, G4double Tmin, G4double Tmax, G4int nbin)
G4PhysicsTableBuildLabTimeTable (G4PhysicsTable *theDEDXTable, G4PhysicsTable *theLabTimeTable, G4double Tmin, G4double Tmax, G4int nbin)
G4PhysicsTableBuildProperTimeTable (G4PhysicsTable *theDEDXTable, G4PhysicsTable *ProperTimeTable, G4double Tmin, G4double Tmax, G4int nbin)
G4PhysicsTableBuildRangeCoeffATable (G4PhysicsTable *theRangeTable, G4PhysicsTable *theCoeffATable, G4double Tmin, G4double Tmax, G4int nbin)
G4PhysicsTableBuildRangeCoeffBTable (G4PhysicsTable *theRangeTable, G4PhysicsTable *theCoeffBTable, G4double Tmin, G4double Tmax, G4int nbin)
G4PhysicsTableBuildRangeCoeffCTable (G4PhysicsTable *theRangeTable, G4PhysicsTable *theCoeffCTable, G4double Tmin, G4double Tmax, G4int nbin)
G4PhysicsTableBuildInverseRangeTable (G4PhysicsTable *theRangeTable, G4PhysicsTable *theRangeCoeffATable, G4PhysicsTable *theRangeCoeffBTable, G4PhysicsTable *theRangeCoeffCTable, G4PhysicsTable *theInverseRangeTable, G4double Tmin, G4double Tmax, G4int nbin)
G4bool CutsWhereModified ()

Static Protected Member Functions

static G4bool EqualCutVectors (G4double *vec1, G4double *vec2)
static G4doubleCopyCutVectors (G4double *dest, G4double *source)

Protected Attributes

G4double ParticleMass

Static Protected Attributes

static G4double dRoverRange = 20*perCent
static G4double finalRange = 1*mm
static G4double finalRangeRequested = -1*mm
static G4double c1lim = dRoverRange
static G4double c2lim = 2.*(1.-dRoverRange)*finalRange
static G4double c3lim = -(1.-dRoverRange)*finalRange*finalRange
static G4bool rndmStepFlag = false
static G4bool EnlossFlucFlag = true
static G4bool subSecFlag = false
static G4double MinDeltaCutInRange = 0.100*mm
static G4doubleMinDeltaEnergy = 0
static G4boolLowerLimitForced = 0
static G4bool setMinDeltaCutInRange = false
static G4EnergyLossMessengerELossMessenger = 0

Detailed Description

Definition at line 68 of file G4VEnergyLoss.hh.


Constructor & Destructor Documentation

G4VEnergyLoss::G4VEnergyLoss ( const G4String ,
G4ProcessType  aType = fNotDefined 
)

Definition at line 81 of file G4VEnergyLoss.cc.

References ELossMessenger, G4cout, G4endl, and ParticleMass.

00082                   : G4VContinuousDiscreteProcess(aName, aType),
00083      lastMaterial(NULL),
00084      nmaxCont1(4),
00085      nmaxCont2(16)
00086 {
00087   if(!ELossMessenger) { 
00088     G4cout << "### G4VEnergyLoss class is obsolete "
00089            << "and will be removed for the next release." << G4endl;
00090     ELossMessenger = new G4EnergyLossMessenger(); 
00091   }
00092 
00093   imat = 0;
00094   f1Fluct = f2Fluct = e1Fluct = e2Fluct = rateFluct = ipotFluct = e1LogFluct 
00095     = e2LogFluct = ipotLogFluct = taulow = tauhigh = ltaulow = ltauhigh 
00096     = ParticleMass = 0.0;
00097 }

G4VEnergyLoss::~G4VEnergyLoss (  )  [virtual]

Definition at line 101 of file G4VEnergyLoss.cc.

00102 {}


Member Function Documentation

virtual G4VParticleChange* G4VEnergyLoss::AlongStepDoIt ( const G4Track track,
const G4Step Step 
) [pure virtual]

Reimplemented from G4VContinuousDiscreteProcess.

G4PhysicsTable * G4VEnergyLoss::BuildInverseRangeTable ( G4PhysicsTable theRangeTable,
G4PhysicsTable theRangeCoeffATable,
G4PhysicsTable theRangeCoeffBTable,
G4PhysicsTable theRangeCoeffCTable,
G4PhysicsTable theInverseRangeTable,
G4double  Tmin,
G4double  Tmax,
G4int  nbin 
) [protected]

Definition at line 618 of file G4VEnergyLoss.cc.

References G4PhysicsTable::clearAndDestroy(), G4PhysicsTable::insert(), and G4PhysicsTable::length().

00626 {
00627   G4double SmallestRange,BiggestRange ;
00628   G4bool isOut ;
00629   size_t numOfMaterials = theRangeTable->length();
00630 
00631     if(theInverseRangeTable)
00632     { theInverseRangeTable->clearAndDestroy();
00633       delete theInverseRangeTable; }
00634     theInverseRangeTable = new G4PhysicsTable(numOfMaterials);
00635 
00636   // loop for materials
00637   for (size_t J=0;  J<numOfMaterials; J++)
00638   {
00639     SmallestRange = (*theRangeTable)(J)->
00640                        GetValue(LowestKineticEnergy,isOut) ;
00641     BiggestRange = (*theRangeTable)(J)->
00642                        GetValue(HighestKineticEnergy,isOut) ;
00643     G4PhysicsLogVector* aVector;
00644     aVector = new G4PhysicsLogVector(SmallestRange,
00645                             BiggestRange,TotBin);
00646 
00647     InvertRangeVector(theRangeTable,
00648                       theRangeCoeffATable,
00649                       theRangeCoeffBTable,
00650                       theRangeCoeffCTable,
00651          LowestKineticEnergy,HighestKineticEnergy,TotBin,J, aVector);
00652 
00653     theInverseRangeTable->insert(aVector);
00654   }
00655   return theInverseRangeTable ;
00656 }

G4PhysicsTable * G4VEnergyLoss::BuildLabTimeTable ( G4PhysicsTable theDEDXTable,
G4PhysicsTable theLabTimeTable,
G4double  Tmin,
G4double  Tmax,
G4int  nbin 
) [protected]

Definition at line 367 of file G4VEnergyLoss.cc.

References G4PhysicsTable::clearAndDestroy(), G4PhysicsTable::insert(), and G4PhysicsTable::length().

00372 {
00373   size_t numOfMaterials = theDEDXTable->length();
00374 
00375   if(theLabTimeTable)
00376   { theLabTimeTable->clearAndDestroy();
00377     delete theLabTimeTable; }
00378   theLabTimeTable = new G4PhysicsTable(numOfMaterials);
00379 
00380 
00381   for (size_t J=0;  J<numOfMaterials; J++)
00382   {
00383     G4PhysicsLogVector* aVector;
00384 
00385     aVector = new G4PhysicsLogVector(LowestKineticEnergy,
00386                             HighestKineticEnergy,TotBin);
00387 
00388     BuildLabTimeVector(theDEDXTable,
00389               LowestKineticEnergy,HighestKineticEnergy,TotBin,J,aVector);
00390     theLabTimeTable->insert(aVector);
00391 
00392 
00393   }
00394   return theLabTimeTable ;
00395 }

G4PhysicsTable * G4VEnergyLoss::BuildProperTimeTable ( G4PhysicsTable theDEDXTable,
G4PhysicsTable ProperTimeTable,
G4double  Tmin,
G4double  Tmax,
G4int  nbin 
) [protected]

Definition at line 399 of file G4VEnergyLoss.cc.

References G4PhysicsTable::clearAndDestroy(), G4PhysicsTable::insert(), and G4PhysicsTable::length().

00404 {
00405   size_t numOfMaterials = theDEDXTable->length();
00406 
00407   if(theProperTimeTable)
00408   { theProperTimeTable->clearAndDestroy();
00409     delete theProperTimeTable; }
00410   theProperTimeTable = new G4PhysicsTable(numOfMaterials);
00411 
00412 
00413   for (size_t J=0;  J<numOfMaterials; J++)
00414   {
00415     G4PhysicsLogVector* aVector;
00416 
00417     aVector = new G4PhysicsLogVector(LowestKineticEnergy,
00418                             HighestKineticEnergy,TotBin);
00419 
00420     BuildProperTimeVector(theDEDXTable,
00421               LowestKineticEnergy,HighestKineticEnergy,TotBin,J,aVector);
00422     theProperTimeTable->insert(aVector);
00423 
00424 
00425   }
00426   return theProperTimeTable ;
00427 }

G4PhysicsTable * G4VEnergyLoss::BuildRangeCoeffATable ( G4PhysicsTable theRangeTable,
G4PhysicsTable theCoeffATable,
G4double  Tmin,
G4double  Tmax,
G4int  nbin 
) [protected]

Definition at line 716 of file G4VEnergyLoss.cc.

References G4PhysicsTable::clearAndDestroy(), G4PhysicsVector::GetValue(), G4PhysicsTable::insert(), G4PhysicsTable::length(), and G4PhysicsVector::PutValue().

00722 {
00723   G4int numOfMaterials = theRangeTable->length();
00724 
00725   if(theRangeCoeffATable)
00726   { theRangeCoeffATable->clearAndDestroy();
00727     delete theRangeCoeffATable; }
00728   theRangeCoeffATable = new G4PhysicsTable(numOfMaterials);
00729 
00730   G4double RTable = std::exp(std::log(HighestKineticEnergy/LowestKineticEnergy)/TotBin) ;
00731   G4double R2 = RTable*RTable ;
00732   G4double R1 = RTable+1.;
00733   G4double w = R1*(RTable-1.)*(RTable-1.);
00734   G4double w1 = RTable/w , w2 = -RTable*R1/w , w3 = R2/w ;
00735   G4double Ti , Tim , Tip , Ri , Rim , Rip , Value ;
00736   G4bool isOut;
00737 
00738   //  loop for materials
00739   for (G4int J=0; J<numOfMaterials; J++)
00740   {
00741     G4int binmax=TotBin ;
00742     G4PhysicsLinearVector* aVector =
00743                            new G4PhysicsLinearVector(0.,binmax, TotBin);
00744     Ti = LowestKineticEnergy ;
00745     G4PhysicsVector* rangeVector= (*theRangeTable)[J];
00746 
00747     for ( G4int i=0; i<TotBin; i++)
00748     {
00749       Ri = rangeVector->GetValue(Ti,isOut) ;
00750       if ( i==0 )
00751         Rim = 0. ;
00752       else
00753       {
00754         Tim = Ti/RTable ;
00755         Rim = rangeVector->GetValue(Tim,isOut);
00756       }
00757       if ( i==(TotBin-1))
00758         Rip = Ri ;
00759       else
00760       {
00761         Tip = Ti*RTable ;
00762         Rip = rangeVector->GetValue(Tip,isOut);
00763       }
00764       Value = (w1*Rip + w2*Ri + w3*Rim)/(Ti*Ti) ;
00765 
00766       aVector->PutValue(i,Value);
00767       Ti = RTable*Ti ;
00768     }
00769  
00770     theRangeCoeffATable->insert(aVector);
00771   }
00772   return theRangeCoeffATable ;
00773 }

G4PhysicsTable * G4VEnergyLoss::BuildRangeCoeffBTable ( G4PhysicsTable theRangeTable,
G4PhysicsTable theCoeffBTable,
G4double  Tmin,
G4double  Tmax,
G4int  nbin 
) [protected]

Definition at line 777 of file G4VEnergyLoss.cc.

References G4PhysicsTable::clearAndDestroy(), G4PhysicsVector::GetValue(), G4PhysicsTable::insert(), G4PhysicsTable::length(), and G4PhysicsVector::PutValue().

00783 {
00784   G4int numOfMaterials = theRangeTable->length();
00785 
00786   if(theRangeCoeffBTable)
00787   { theRangeCoeffBTable->clearAndDestroy();
00788     delete theRangeCoeffBTable; }
00789   theRangeCoeffBTable = new G4PhysicsTable(numOfMaterials);
00790 
00791   G4double RTable = std::exp(std::log(HighestKineticEnergy/LowestKineticEnergy)/TotBin) ;
00792   G4double R2 = RTable*RTable ;
00793   G4double R1 = RTable+1.;
00794   G4double w = R1*(RTable-1.)*(RTable-1.);
00795   G4double w1 = -R1/w , w2 = R1*(R2+1.)/w , w3 = -R2*R1/w ;
00796   G4double Ti , Tim , Tip , Ri , Rim , Rip , Value ;
00797   G4bool isOut;
00798 
00799   //  loop for materials
00800   for (G4int J=0; J<numOfMaterials; J++)
00801   {
00802     G4int binmax=TotBin ;
00803     G4PhysicsLinearVector* aVector =
00804                         new G4PhysicsLinearVector(0.,binmax, TotBin);
00805     Ti = LowestKineticEnergy ;
00806     G4PhysicsVector* rangeVector= (*theRangeTable)[J];
00807   
00808     for ( G4int i=0; i<TotBin; i++)
00809     {
00810       Ri = rangeVector->GetValue(Ti,isOut) ;
00811       if ( i==0 )
00812          Rim = 0. ;
00813       else
00814       {
00815         Tim = Ti/RTable ;
00816         Rim = rangeVector->GetValue(Tim,isOut);
00817       }
00818       if ( i==(TotBin-1))
00819         Rip = Ri ;
00820       else
00821       {
00822         Tip = Ti*RTable ;
00823         Rip = rangeVector->GetValue(Tip,isOut);
00824       }
00825       Value = (w1*Rip + w2*Ri + w3*Rim)/Ti;
00826 
00827       aVector->PutValue(i,Value);
00828       Ti = RTable*Ti ;
00829     }
00830     theRangeCoeffBTable->insert(aVector);
00831   }
00832   return theRangeCoeffBTable ;
00833 }

G4PhysicsTable * G4VEnergyLoss::BuildRangeCoeffCTable ( G4PhysicsTable theRangeTable,
G4PhysicsTable theCoeffCTable,
G4double  Tmin,
G4double  Tmax,
G4int  nbin 
) [protected]

Definition at line 837 of file G4VEnergyLoss.cc.

References G4PhysicsTable::clearAndDestroy(), G4PhysicsVector::GetValue(), G4PhysicsTable::insert(), G4PhysicsTable::length(), and G4PhysicsVector::PutValue().

00843 {
00844   G4int numOfMaterials = theRangeTable->length();
00845 
00846   if(theRangeCoeffCTable)
00847   { theRangeCoeffCTable->clearAndDestroy();
00848     delete theRangeCoeffCTable; }
00849   theRangeCoeffCTable = new G4PhysicsTable(numOfMaterials);
00850 
00851   G4double RTable = std::exp(std::log(HighestKineticEnergy/LowestKineticEnergy)/TotBin) ;
00852   G4double R2 = RTable*RTable ;
00853   G4double R1 = RTable+1.;
00854   G4double w = R1*(RTable-1.)*(RTable-1.);
00855   G4double w1 = 1./w , w2 = -RTable*R1/w , w3 = RTable*R2/w ;
00856   G4double Ti , Tim , Tip , Ri , Rim , Rip , Value ;
00857   G4bool isOut;
00858 
00859   //  loop for materials
00860   for (G4int J=0; J<numOfMaterials; J++)
00861   {
00862     G4int binmax=TotBin ;
00863     G4PhysicsLinearVector* aVector =
00864                       new G4PhysicsLinearVector(0.,binmax, TotBin);
00865     Ti = LowestKineticEnergy ;
00866     G4PhysicsVector* rangeVector= (*theRangeTable)[J];
00867   
00868     for ( G4int i=0; i<TotBin; i++)
00869     {
00870       Ri = rangeVector->GetValue(Ti,isOut) ;
00871       if ( i==0 )
00872         Rim = 0. ;
00873       else
00874       {
00875         Tim = Ti/RTable ;
00876         Rim = rangeVector->GetValue(Tim,isOut);
00877       }
00878       if ( i==(TotBin-1))
00879         Rip = Ri ;
00880       else
00881       {
00882         Tip = Ti*RTable ;
00883         Rip = rangeVector->GetValue(Tip,isOut);
00884       }
00885       Value = w1*Rip + w2*Ri + w3*Rim ;
00886 
00887       aVector->PutValue(i,Value);
00888       Ti = RTable*Ti ;
00889     }
00890     theRangeCoeffCTable->insert(aVector);
00891   }
00892   return theRangeCoeffCTable ;
00893 }

G4PhysicsTable * G4VEnergyLoss::BuildRangeTable ( G4PhysicsTable theDEDXTable,
G4PhysicsTable theRangeTable,
G4double  Tmin,
G4double  Tmax,
G4int  nbin 
) [protected]

Definition at line 133 of file G4VEnergyLoss.cc.

References G4PhysicsTable::clearAndDestroy(), G4PhysicsTable::insert(), and G4PhysicsTable::length().

00137 {
00138    size_t numOfMaterials = theDEDXTable->length();
00139 
00140    if(theRangeTable)
00141    { theRangeTable->clearAndDestroy();
00142      delete theRangeTable; }
00143    theRangeTable = new G4PhysicsTable(numOfMaterials);
00144 
00145    // loop for materials
00146 
00147    for (size_t J=0;  J<numOfMaterials; J++)
00148    {
00149      G4PhysicsLogVector* aVector;
00150      aVector = new G4PhysicsLogVector(LowestKineticEnergy,
00151                               HighestKineticEnergy,TotBin);
00152 
00153      BuildRangeVector(theDEDXTable,LowestKineticEnergy,HighestKineticEnergy,
00154                       TotBin,J,aVector);
00155      theRangeTable->insert(aVector);
00156    }
00157    return theRangeTable ;
00158 }

G4double * G4VEnergyLoss::CopyCutVectors ( G4double dest,
G4double source 
) [static, protected]

Definition at line 1131 of file G4VEnergyLoss.cc.

References G4Material::GetNumberOfMaterials().

01132 {
01133   if ( dest != 0) delete [] dest;
01134   dest = new G4double [G4Material::GetNumberOfMaterials()];
01135   for (size_t j=0; j<G4Material::GetNumberOfMaterials(); j++){
01136     dest[j] = source[j];
01137   }
01138   return dest;
01139 }

G4bool G4VEnergyLoss::CutsWhereModified (  )  [protected]

Definition at line 1143 of file G4VEnergyLoss.cc.

References G4ProductionCutsTable::GetMaterialCutsCouple(), G4ProductionCutsTable::GetProductionCutsTable(), G4ProductionCutsTable::GetTableSize(), and G4MaterialCutsCouple::IsRecalcNeeded().

01144 {
01145   G4bool wasModified = false;
01146   const G4ProductionCutsTable* theCoupleTable=
01147         G4ProductionCutsTable::GetProductionCutsTable();
01148   size_t numOfCouples = theCoupleTable->GetTableSize();
01149 
01150   for (size_t j=0; j<numOfCouples; j++){
01151     if (theCoupleTable->GetMaterialCutsCouple(j)->IsRecalcNeeded()) {
01152       wasModified = true;
01153       break;
01154     }
01155   }
01156   return wasModified;
01157 }

G4bool G4VEnergyLoss::EqualCutVectors ( G4double vec1,
G4double vec2 
) [static, protected]

Definition at line 1116 of file G4VEnergyLoss.cc.

References G4Material::GetNumberOfMaterials().

01117 {
01118   if ( (vec1==0 ) || (vec2==0) ) return false;
01119   
01120   G4bool flag = true;
01121    
01122   for (size_t j=0; flag && j<G4Material::GetNumberOfMaterials(); j++){
01123     flag = (vec1[j] == vec2[j]);
01124   }
01125 
01126   return flag;
01127 }

virtual G4double G4VEnergyLoss::GetContinuousStepLimit ( const G4Track track,
G4double  previousStepSize,
G4double  currentMinimumStep,
G4double currentSafety 
) [pure virtual]

Implements G4VContinuousDiscreteProcess.

G4double G4VEnergyLoss::GetLossWithFluct ( const G4DynamicParticle aParticle,
const G4MaterialCutsCouple couple,
G4double  ChargeSquare,
G4double  MeanLoss,
G4double  step 
) [protected]

Definition at line 897 of file G4VEnergyLoss.cc.

References G4Poisson(), G4UniformRand, G4Material::GetElectronDensity(), G4IonisParamMat::GetEnergy0fluct(), G4IonisParamMat::GetEnergy1fluct(), G4IonisParamMat::GetEnergy2fluct(), G4IonisParamMat::GetF1fluct(), G4IonisParamMat::GetF2fluct(), G4MaterialCutsCouple::GetIndex(), G4Material::GetIonisation(), G4DynamicParticle::GetKineticEnergy(), G4IonisParamMat::GetLogEnergy1fluct(), G4IonisParamMat::GetLogEnergy2fluct(), G4IonisParamMat::GetLogMeanExcEnergy(), G4DynamicParticle::GetMass(), G4MaterialCutsCouple::GetMaterial(), G4IonisParamMat::GetMeanExcitationEnergy(), G4ProductionCutsTable::GetProductionCutsTable(), G4IonisParamMat::GetRateionexcfluct(), and ParticleMass.

00904 {
00905   const G4double minLoss = 1.*eV ;
00906   const G4double probLim = 0.01 ;
00907   const G4double sumaLim = -std::log(probLim) ;
00908   const G4double alim=10.;
00909   const G4double kappa = 10. ;
00910   const G4double factor = twopi_mc2_rcl2 ;
00911   const G4Material* aMaterial = couple->GetMaterial();
00912 
00913   // check if the material has changed ( cache mechanism)
00914 
00915   if (aMaterial != lastMaterial)
00916     {
00917       lastMaterial = aMaterial;
00918       imat         = couple->GetIndex();
00919       f1Fluct      = aMaterial->GetIonisation()->GetF1fluct();
00920       f2Fluct      = aMaterial->GetIonisation()->GetF2fluct();
00921       e1Fluct      = aMaterial->GetIonisation()->GetEnergy1fluct();
00922       e2Fluct      = aMaterial->GetIonisation()->GetEnergy2fluct();
00923       e1LogFluct   = aMaterial->GetIonisation()->GetLogEnergy1fluct();
00924       e2LogFluct   = aMaterial->GetIonisation()->GetLogEnergy2fluct();
00925       rateFluct    = aMaterial->GetIonisation()->GetRateionexcfluct();
00926       ipotFluct    = aMaterial->GetIonisation()->GetMeanExcitationEnergy();
00927       ipotLogFluct = aMaterial->GetIonisation()->GetLogMeanExcEnergy();
00928     }
00929   G4double threshold,w1,w2,C,
00930            beta2,suma,e0,loss,lossc ,w,electronDensity;
00931   G4double a1,a2,a3;
00932   G4double p1,p2,p3 ;
00933   G4int nb;
00934   G4double Corrfac, na,alfa,rfac,namean,sa,alfa1,ea,sea;
00935   G4double dp3;
00936   G4double siga ;
00937 
00938   // shortcut for very very small loss
00939   if(MeanLoss < minLoss) return MeanLoss ;
00940 
00941   // get particle data
00942   G4double Tkin   = aParticle->GetKineticEnergy();
00943   ParticleMass = aParticle->GetMass() ;
00944 
00945   threshold = (*((G4ProductionCutsTable::GetProductionCutsTable())
00946                 ->GetEnergyCutsVector(1)))[imat];
00947   G4double rmass = electron_mass_c2/ParticleMass;
00948   G4double tau   = Tkin/ParticleMass, tau1 = tau+1., tau2 = tau*(tau+2.);
00949   G4double Tm    = 2.*electron_mass_c2*tau2/(1.+2.*tau1*rmass+rmass*rmass);
00950 
00951   if(Tm > threshold) Tm = threshold;
00952 
00953   beta2 = tau2/(tau1*tau1);
00954 
00955   // Gaussian fluctuation ?
00956   if(MeanLoss >= kappa*Tm || Tm < kappa*ipotFluct)
00957   {
00958     electronDensity = aMaterial->GetElectronDensity() ;
00959     siga = std::sqrt(Tm*(1.0-0.5*beta2)*step*
00960                 factor*electronDensity*ChargeSquare/beta2) ;
00961     do {
00962       loss = G4RandGauss::shoot(MeanLoss,siga) ;
00963     } while (loss < 0. || loss > 2.*MeanLoss);
00964     return loss;
00965   }
00966 
00967   w1 = Tm/ipotFluct;
00968   w2 = std::log(2.*electron_mass_c2*tau2);
00969 
00970   C = MeanLoss*(1.-rateFluct)/(w2-ipotLogFluct-beta2);
00971 
00972   a1 = C*f1Fluct*(w2-e1LogFluct-beta2)/e1Fluct;
00973   a2 = C*f2Fluct*(w2-e2LogFluct-beta2)/e2Fluct;
00974   a3 = rateFluct*MeanLoss*(Tm-ipotFluct)/(ipotFluct*Tm*std::log(w1));
00975   if(a1 < 0.) a1 = 0.;
00976   if(a2 < 0.) a2 = 0.;
00977   if(a3 < 0.) a3 = 0.;
00978 
00979   suma = a1+a2+a3;
00980 
00981   loss = 0. ;
00982 
00983   if(suma < sumaLim)             // very small Step
00984     {
00985       e0 = aMaterial->GetIonisation()->GetEnergy0fluct();
00986 
00987       if(Tm == ipotFluct)
00988       {
00989         a3 = MeanLoss/e0;
00990 
00991         if(a3>alim)
00992         {
00993           siga=std::sqrt(a3) ;
00994           p3 = std::max(0.,G4RandGauss::shoot(a3,siga)+0.5);
00995         }
00996         else
00997           p3 = G4float(G4Poisson(a3));
00998 
00999         loss = p3*e0 ;
01000 
01001         if(p3 > 0)
01002           loss += (1.-2.*G4UniformRand())*e0 ;
01003 
01004       }
01005       else
01006       {
01007         Tm = Tm-ipotFluct+e0 ;
01008         a3 = MeanLoss*(Tm-e0)/(Tm*e0*std::log(Tm/e0));
01009 
01010         if(a3>alim)
01011         {
01012           siga=std::sqrt(a3) ;
01013           p3 = std::max(0.,G4RandGauss::shoot(a3,siga)+0.5);
01014         }
01015         else
01016           p3 = G4float(G4Poisson(a3));
01017 
01018         if(p3 > 0)
01019         {
01020           w = (Tm-e0)/Tm ;
01021           if(p3 > G4float(nmaxCont2))
01022           {
01023             dp3 = p3 ;
01024             Corrfac = dp3/G4float(nmaxCont2) ;
01025             p3 = G4float(nmaxCont2) ;
01026           }
01027           else
01028             Corrfac = 1. ;
01029 
01030           for(G4int i=0; i<G4int(p3); i++) loss += 1./(1.-w*G4UniformRand()) ;
01031           loss *= e0*Corrfac ;  
01032         }        
01033       }
01034     }
01035     
01036   else                              // not so small Step
01037     {
01038       // excitation type 1
01039       if(a1>alim)
01040       {
01041         siga=std::sqrt(a1) ;
01042         p1 = std::max(0.,G4RandGauss::shoot(a1,siga)+0.5);
01043       }
01044       else
01045        p1 = G4float(G4Poisson(a1));
01046 
01047       // excitation type 2
01048       if(a2>alim)
01049       {
01050         siga=std::sqrt(a2) ;
01051         p2 = std::max(0.,G4RandGauss::shoot(a2,siga)+0.5);
01052       }
01053       else
01054         p2 = G4float(G4Poisson(a2));
01055 
01056       loss = p1*e1Fluct+p2*e2Fluct;
01057       // smearing to avoid unphysical peaks
01058       if(p2 > 0)
01059         loss += (1.-2.*G4UniformRand())*e2Fluct;
01060       else if (loss>0.)
01061         loss += (1.-2.*G4UniformRand())*e1Fluct;   
01062 
01063       // ionisation .......................................
01064      if(a3 > 0.)
01065      {
01066       if(a3>alim)
01067       {
01068         siga=std::sqrt(a3) ;
01069         p3 = std::max(0.,G4RandGauss::shoot(a3,siga)+0.5);
01070       }
01071       else
01072         p3 = G4float(G4Poisson(a3));
01073 
01074       lossc = 0.;
01075       if(p3 > 0)
01076       {
01077         na = 0.; 
01078         alfa = 1.;
01079         if (p3 > G4float(nmaxCont2))
01080         {
01081           dp3        = p3;
01082           rfac       = dp3/(G4float(nmaxCont2)+dp3);
01083           namean     = p3*rfac;
01084           sa         = G4float(nmaxCont1)*rfac;
01085           na         = G4RandGauss::shoot(namean,sa);
01086           if (na > 0.)
01087           {
01088             alfa   = w1*(G4float(nmaxCont2)+p3)/(w1*G4float(nmaxCont2)+p3);
01089             alfa1  = alfa*std::log(alfa)/(alfa-1.);
01090             ea     = na*ipotFluct*alfa1;
01091             sea    = ipotFluct*std::sqrt(na*(alfa-alfa1*alfa1));
01092             lossc += G4RandGauss::shoot(ea,sea);
01093           }
01094         }
01095 
01096         nb = G4int(p3-na);
01097         if (nb > 0)
01098         {
01099           w2 = alfa*ipotFluct;
01100           w  = (Tm-w2)/Tm;      
01101           for (G4int k=0; k<nb; k++) lossc += w2/(1.-w*G4UniformRand());
01102         }
01103       }        
01104       loss += lossc;  
01105      }
01106     } 
01107 
01108   if( loss < 0.)
01109     loss = 0.;
01110 
01111   return loss ;
01112 }

virtual G4double G4VEnergyLoss::GetMeanFreePath ( const G4Track track,
G4double  previousStepSize,
G4ForceCondition condition 
) [pure virtual]

Implements G4VContinuousDiscreteProcess.

virtual G4VParticleChange* G4VEnergyLoss::PostStepDoIt ( const G4Track track,
const G4Step Step 
) [pure virtual]

Reimplemented from G4VContinuousDiscreteProcess.

void G4VEnergyLoss::SetEnlossFluc ( G4bool  value  )  [static]

Definition at line 110 of file G4VEnergyLoss.cc.

References EnlossFlucFlag.

00110 {EnlossFlucFlag = value;}

void G4VEnergyLoss::SetMinDeltaCutInRange ( G4double  value  )  [static]

Definition at line 118 of file G4VEnergyLoss.cc.

References MinDeltaCutInRange, and setMinDeltaCutInRange.

00119 {MinDeltaCutInRange = value; setMinDeltaCutInRange = true;}

void G4VEnergyLoss::SetRndmStep ( G4bool  value  )  [static]

Definition at line 106 of file G4VEnergyLoss.cc.

References rndmStepFlag.

00106 {rndmStepFlag = value;}

void G4VEnergyLoss::SetStepFunction ( G4double  c1,
G4double  c2 
) [static]

Definition at line 123 of file G4VEnergyLoss.cc.

References c1lim, c2lim, c3lim, dRoverRange, and finalRangeRequested.

void G4VEnergyLoss::SetSubSec ( G4bool  value  )  [static]

Definition at line 114 of file G4VEnergyLoss.cc.

References subSecFlag.

00114 {subSecFlag = value;}


Field Documentation

G4double G4VEnergyLoss::c1lim = dRoverRange [static, protected]

Definition at line 231 of file G4VEnergyLoss.hh.

Referenced by SetStepFunction().

G4double G4VEnergyLoss::c2lim = 2.*(1.-dRoverRange)*finalRange [static, protected]

Definition at line 231 of file G4VEnergyLoss.hh.

Referenced by SetStepFunction().

G4double G4VEnergyLoss::c3lim = -(1.-dRoverRange)*finalRange*finalRange [static, protected]

Definition at line 231 of file G4VEnergyLoss.hh.

Referenced by SetStepFunction().

G4double G4VEnergyLoss::dRoverRange = 20*perCent [static, protected]

Definition at line 227 of file G4VEnergyLoss.hh.

Referenced by SetStepFunction().

G4EnergyLossMessenger * G4VEnergyLoss::ELossMessenger = 0 [static, protected]

Definition at line 243 of file G4VEnergyLoss.hh.

Referenced by G4VEnergyLoss().

G4bool G4VEnergyLoss::EnlossFlucFlag = true [static, protected]

Definition at line 234 of file G4VEnergyLoss.hh.

Referenced by SetEnlossFluc().

G4double G4VEnergyLoss::finalRange = 1*mm [static, protected]

Definition at line 229 of file G4VEnergyLoss.hh.

G4double G4VEnergyLoss::finalRangeRequested = -1*mm [static, protected]

Definition at line 230 of file G4VEnergyLoss.hh.

Referenced by SetStepFunction().

G4bool * G4VEnergyLoss::LowerLimitForced = 0 [static, protected]

Definition at line 239 of file G4VEnergyLoss.hh.

G4double G4VEnergyLoss::MinDeltaCutInRange = 0.100*mm [static, protected]

Definition at line 237 of file G4VEnergyLoss.hh.

Referenced by SetMinDeltaCutInRange().

G4double * G4VEnergyLoss::MinDeltaEnergy = 0 [static, protected]

Definition at line 238 of file G4VEnergyLoss.hh.

G4double G4VEnergyLoss::ParticleMass [protected]

Definition at line 175 of file G4VEnergyLoss.hh.

Referenced by G4VEnergyLoss(), and GetLossWithFluct().

G4bool G4VEnergyLoss::rndmStepFlag = false [static, protected]

Definition at line 233 of file G4VEnergyLoss.hh.

Referenced by SetRndmStep().

G4bool G4VEnergyLoss::setMinDeltaCutInRange = false [static, protected]

Definition at line 241 of file G4VEnergyLoss.hh.

Referenced by SetMinDeltaCutInRange().

G4bool G4VEnergyLoss::subSecFlag = false [static, protected]

Definition at line 235 of file G4VEnergyLoss.hh.

Referenced by SetSubSec().


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