G4hRDEnergyLoss Class Reference

#include <G4hRDEnergyLoss.hh>

Inheritance diagram for G4hRDEnergyLoss:

G4VContinuousDiscreteProcess G4VProcess G4hImpactIonisation

Public Member Functions

 G4hRDEnergyLoss (const G4String &)
 ~G4hRDEnergyLoss ()
virtual G4double GetMeanFreePath (const G4Track &track, G4double previousStepSize, enum G4ForceCondition *condition)=0
virtual G4VParticleChangePostStepDoIt (const G4Track &track, const G4Step &Step)=0

Static Public Member Functions

static G4int GetNumberOfProcesses ()
static void SetNumberOfProcesses (G4int number)
static void PlusNumberOfProcesses ()
static void MinusNumberOfProcesses ()
static void SetdRoverRange (G4double value)
static void SetRndmStep (G4bool value)
static void SetEnlossFluc (G4bool value)
static void SetStepFunction (G4double c1, G4double c2)

Protected Member Functions

G4bool CutsWhereModified ()

Static Protected Member Functions

static void BuildDEDXTable (const G4ParticleDefinition &aParticleType)

Protected Attributes

const G4double MaxExcitationNumber
const G4double probLimFluct
const long nmaxDirectFluct
const long nmaxCont1
const long nmaxCont2
G4PhysicsTabletheLossTable
G4double linLossLimit
G4double MinKineticEnergy

Static Protected Attributes

static G4PhysicsTabletheDEDXpTable = 0
static G4PhysicsTabletheDEDXpbarTable = 0
static G4PhysicsTabletheRangepTable = 0
static G4PhysicsTabletheRangepbarTable = 0
static G4PhysicsTabletheInverseRangepTable = 0
static G4PhysicsTabletheInverseRangepbarTable = 0
static G4PhysicsTabletheLabTimepTable = 0
static G4PhysicsTabletheLabTimepbarTable = 0
static G4PhysicsTabletheProperTimepTable = 0
static G4PhysicsTabletheProperTimepbarTable = 0
static G4PhysicsTable ** RecorderOfpProcess
static G4PhysicsTable ** RecorderOfpbarProcess
static G4int CounterOfpProcess = 0
static G4int CounterOfpbarProcess = 0
static G4double ParticleMass
static G4double ptableElectronCutInRange = 0.0*mm
static G4double pbartableElectronCutInRange = 0.0*mm
static G4double Charge
static G4double LowestKineticEnergy = 10.*eV
static G4double HighestKineticEnergy = 100.*GeV
static G4int TotBin = 360
static G4double RTable
static G4double LOGRTable
static G4double dRoverRange = 0.20
static G4double finalRange = 200.*micrometer
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

Detailed Description

Definition at line 93 of file G4hRDEnergyLoss.hh.


Constructor & Destructor Documentation

G4hRDEnergyLoss::G4hRDEnergyLoss ( const G4String  ) 

Definition at line 163 of file G4hRDEnergyLoss.cc.

00164   : G4VContinuousDiscreteProcess (processName),
00165     MaxExcitationNumber (1.e6),
00166     probLimFluct (0.01),
00167     nmaxDirectFluct (100),
00168     nmaxCont1(4),
00169     nmaxCont2(16),
00170     theLossTable(0),
00171     linLossLimit(0.05),
00172     MinKineticEnergy(0.0) 
00173 {;}

G4hRDEnergyLoss::~G4hRDEnergyLoss (  ) 

Definition at line 177 of file G4hRDEnergyLoss.cc.

References G4PhysicsTable::clearAndDestroy(), and theLossTable.

00178 {
00179   if(theLossTable) {
00180     theLossTable->clearAndDestroy();
00181     delete theLossTable;
00182   }
00183 }


Member Function Documentation

void G4hRDEnergyLoss::BuildDEDXTable ( const G4ParticleDefinition aParticleType  )  [static, protected]

Definition at line 247 of file G4hRDEnergyLoss.cc.

References Charge, CounterOfpbarProcess, CounterOfpProcess, G4ParticleDefinition::GetPDGCharge(), G4ParticleDefinition::GetPDGMass(), G4ProductionCutsTable::GetProductionCutsTable(), G4ProductionCutsTable::GetTableSize(), HighestKineticEnergy, G4PhysicsTable::insert(), LowestKineticEnergy, ParticleMass, RecorderOfpbarProcess, RecorderOfpProcess, G4EnergyLossTables::Register(), theDEDXpbarTable, theDEDXpTable, theInverseRangepbarTable, theInverseRangepTable, theLabTimepbarTable, theLabTimepTable, theProperTimepbarTable, theProperTimepTable, theRangepbarTable, theRangepTable, and TotBin.

Referenced by G4hImpactIonisation::BuildPhysicsTable().

00249 {
00250   //  calculate data members TotBin,LOGRTable,RTable first
00251 
00252   const G4ProductionCutsTable* theCoupleTable=
00253     G4ProductionCutsTable::GetProductionCutsTable();
00254   size_t numOfCouples = theCoupleTable->GetTableSize();
00255 
00256   // create/fill proton or antiproton tables depending on the charge
00257   Charge = aParticleType.GetPDGCharge()/eplus;
00258   ParticleMass = aParticleType.GetPDGMass() ;
00259 
00260   if (Charge>0.) {theDEDXTable= theDEDXpTable;}
00261   else           {theDEDXTable= theDEDXpbarTable;}
00262 
00263   if( ((Charge>0.) && (theDEDXTable==0)) ||
00264       ((Charge<0.) && (theDEDXTable==0)) 
00265       )
00266     {
00267 
00268       // Build energy loss table as a sum of the energy loss due to the
00269       //              different processes.
00270       if( Charge >0.)
00271         {
00272           RecorderOfProcess=RecorderOfpProcess;
00273           CounterOfProcess=CounterOfpProcess;
00274 
00275           if(CounterOfProcess == NumberOfProcesses)
00276             {
00277               if(theDEDXpTable)
00278                 { theDEDXpTable->clearAndDestroy();
00279                   delete theDEDXpTable; }
00280               theDEDXpTable = new G4PhysicsTable(numOfCouples);
00281               theDEDXTable = theDEDXpTable;
00282             }
00283         }
00284       else
00285         {
00286           RecorderOfProcess=RecorderOfpbarProcess;
00287           CounterOfProcess=CounterOfpbarProcess;
00288 
00289           if(CounterOfProcess == NumberOfProcesses)
00290             {
00291               if(theDEDXpbarTable)
00292                 { theDEDXpbarTable->clearAndDestroy();
00293                   delete theDEDXpbarTable; }
00294               theDEDXpbarTable = new G4PhysicsTable(numOfCouples);
00295               theDEDXTable = theDEDXpbarTable;
00296             }
00297         }
00298 
00299       if(CounterOfProcess == NumberOfProcesses)
00300         {
00301           //  loop for materials
00302           G4double LowEdgeEnergy , Value ;
00303           G4bool isOutRange ;
00304           G4PhysicsTable* pointer ;
00305 
00306           for (size_t J=0; J<numOfCouples; J++)
00307             {
00308               // create physics vector and fill it
00309               G4PhysicsLogVector* aVector = new G4PhysicsLogVector(
00310                                                                    LowestKineticEnergy, HighestKineticEnergy, TotBin);
00311 
00312               // loop for the kinetic energy
00313               for (G4int i=0; i<TotBin; i++)
00314                 {
00315                   LowEdgeEnergy = aVector->GetLowEdgeEnergy(i) ;
00316                   Value = 0. ;
00317 
00318                   // loop for the contributing processes
00319                   for (G4int process=0; process < NumberOfProcesses; process++)
00320                     {
00321                       pointer= RecorderOfProcess[process];
00322                       Value += (*pointer)[J]->
00323                         GetValue(LowEdgeEnergy,isOutRange) ;
00324                     }
00325 
00326                   aVector->PutValue(i,Value) ;
00327                 }
00328 
00329               theDEDXTable->insert(aVector) ;
00330             }
00331 
00332           //  reset counter to zero ..................
00333           if( Charge >0.)
00334             CounterOfpProcess=0 ;
00335           else
00336             CounterOfpbarProcess=0 ;
00337 
00338           // Build range table
00339           BuildRangeTable( aParticleType);
00340 
00341           // Build lab/proper time tables
00342           BuildTimeTables( aParticleType) ;
00343 
00344           // Build coeff tables for the energy loss calculation
00345           BuildRangeCoeffATable( aParticleType);
00346           BuildRangeCoeffBTable( aParticleType);
00347           BuildRangeCoeffCTable( aParticleType);
00348 
00349           // invert the range table
00350 
00351           BuildInverseRangeTable(aParticleType);
00352         }
00353     }
00354   // make the energy loss and the range table available
00355 
00356   G4EnergyLossTables::Register(&aParticleType,
00357                                (Charge>0)?
00358                                theDEDXpTable: theDEDXpbarTable,
00359                                (Charge>0)?
00360                                theRangepTable: theRangepbarTable,
00361                                (Charge>0)?
00362                                theInverseRangepTable: theInverseRangepbarTable,
00363                                (Charge>0)?
00364                                theLabTimepTable: theLabTimepbarTable,
00365                                (Charge>0)?
00366                                theProperTimepTable: theProperTimepbarTable,
00367                                LowestKineticEnergy, HighestKineticEnergy,
00368                                proton_mass_c2/aParticleType.GetPDGMass(),
00369                                TotBin);
00370 
00371 }

G4bool G4hRDEnergyLoss::CutsWhereModified (  )  [protected]

Definition at line 1141 of file G4hRDEnergyLoss.cc.

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

Referenced by G4hImpactIonisation::BuildPhysicsTable().

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

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

Implements G4VContinuousDiscreteProcess.

Implemented in G4hImpactIonisation.

G4int G4hRDEnergyLoss::GetNumberOfProcesses (  )  [static]

Definition at line 187 of file G4hRDEnergyLoss.cc.

00188 {   
00189   return NumberOfProcesses; 
00190 } 

void G4hRDEnergyLoss::MinusNumberOfProcesses (  )  [static]

Definition at line 208 of file G4hRDEnergyLoss.cc.

00209 { 
00210   NumberOfProcesses--; 
00211 } 

void G4hRDEnergyLoss::PlusNumberOfProcesses (  )  [static]

Definition at line 201 of file G4hRDEnergyLoss.cc.

00202 { 
00203   NumberOfProcesses++; 
00204 } 

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

Reimplemented from G4VContinuousDiscreteProcess.

Implemented in G4hImpactIonisation.

void G4hRDEnergyLoss::SetdRoverRange ( G4double  value  )  [static]

Definition at line 215 of file G4hRDEnergyLoss.cc.

References dRoverRange.

00216 {
00217   dRoverRange = value;
00218 }

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

Definition at line 229 of file G4hRDEnergyLoss.cc.

References EnlossFlucFlag.

00230 {
00231   EnlossFlucFlag = value;
00232 }

void G4hRDEnergyLoss::SetNumberOfProcesses ( G4int  number  )  [static]

Definition at line 194 of file G4hRDEnergyLoss.cc.

00195 {
00196   NumberOfProcesses=number; 
00197 } 

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

Definition at line 222 of file G4hRDEnergyLoss.cc.

References rndmStepFlag.

00223 {
00224   rndmStepFlag = value;
00225 }

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

Definition at line 236 of file G4hRDEnergyLoss.cc.

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

00237 {
00238   dRoverRange = c1; 
00239   finalRange = c2;
00240   c1lim=dRoverRange;
00241   c2lim=2.*(1-dRoverRange)*finalRange;
00242   c3lim=-(1.-dRoverRange)*finalRange*finalRange;
00243 }


Field Documentation

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

Definition at line 194 of file G4hRDEnergyLoss.hh.

Referenced by SetStepFunction().

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

Definition at line 194 of file G4hRDEnergyLoss.hh.

Referenced by SetStepFunction().

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

Definition at line 194 of file G4hRDEnergyLoss.hh.

Referenced by SetStepFunction().

G4double G4hRDEnergyLoss::Charge [static, protected]

Definition at line 174 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable().

G4int G4hRDEnergyLoss::CounterOfpbarProcess = 0 [static, protected]

Definition at line 165 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable(), and G4hImpactIonisation::BuildPhysicsTable().

G4int G4hRDEnergyLoss::CounterOfpProcess = 0 [static, protected]

Definition at line 164 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable(), and G4hImpactIonisation::BuildPhysicsTable().

G4double G4hRDEnergyLoss::dRoverRange = 0.20 [static, protected]

Definition at line 191 of file G4hRDEnergyLoss.hh.

Referenced by SetdRoverRange(), and SetStepFunction().

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

Definition at line 197 of file G4hRDEnergyLoss.hh.

Referenced by G4hImpactIonisation::AlongStepDoIt(), and SetEnlossFluc().

G4double G4hRDEnergyLoss::finalRange = 200.*micrometer [static, protected]

Definition at line 193 of file G4hRDEnergyLoss.hh.

Referenced by SetStepFunction().

G4double G4hRDEnergyLoss::HighestKineticEnergy = 100.*GeV [static, protected]

Definition at line 177 of file G4hRDEnergyLoss.hh.

Referenced by G4hImpactIonisation::AlongStepDoIt(), BuildDEDXTable(), G4hImpactIonisation::BuildPhysicsTable(), G4hImpactIonisation::GetMeanFreePath(), and G4hImpactIonisation::PrintInfoDefinition().

G4double G4hRDEnergyLoss::linLossLimit [protected]

Definition at line 187 of file G4hRDEnergyLoss.hh.

Referenced by G4hImpactIonisation::AlongStepDoIt().

G4double G4hRDEnergyLoss::LOGRTable [static, protected]

Definition at line 181 of file G4hRDEnergyLoss.hh.

G4double G4hRDEnergyLoss::LowestKineticEnergy = 10.*eV [static, protected]

Definition at line 176 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable(), and G4hImpactIonisation::BuildPhysicsTable().

const G4double G4hRDEnergyLoss::MaxExcitationNumber [protected]

Definition at line 136 of file G4hRDEnergyLoss.hh.

G4double G4hRDEnergyLoss::MinKineticEnergy [protected]

Definition at line 189 of file G4hRDEnergyLoss.hh.

Referenced by G4hImpactIonisation::AlongStepDoIt(), and G4hImpactIonisation::PostStepDoIt().

const long G4hRDEnergyLoss::nmaxCont1 [protected]

Definition at line 138 of file G4hRDEnergyLoss.hh.

const long G4hRDEnergyLoss::nmaxCont2 [protected]

Definition at line 138 of file G4hRDEnergyLoss.hh.

const long G4hRDEnergyLoss::nmaxDirectFluct [protected]

Definition at line 138 of file G4hRDEnergyLoss.hh.

G4double G4hRDEnergyLoss::ParticleMass [static, protected]

Definition at line 168 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable().

G4double G4hRDEnergyLoss::pbartableElectronCutInRange = 0.0*mm [static, protected]

Definition at line 172 of file G4hRDEnergyLoss.hh.

const G4double G4hRDEnergyLoss::probLimFluct [protected]

Definition at line 137 of file G4hRDEnergyLoss.hh.

G4double G4hRDEnergyLoss::ptableElectronCutInRange = 0.0*mm [static, protected]

Definition at line 171 of file G4hRDEnergyLoss.hh.

G4PhysicsTable ** G4hRDEnergyLoss::RecorderOfpbarProcess [static, protected]

Initial value:

  new G4PhysicsTable*[100]

Definition at line 163 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable(), and G4hImpactIonisation::BuildPhysicsTable().

G4PhysicsTable ** G4hRDEnergyLoss::RecorderOfpProcess [static, protected]

Initial value:

  new G4PhysicsTable*[100]

Definition at line 162 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable(), and G4hImpactIonisation::BuildPhysicsTable().

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

Definition at line 196 of file G4hRDEnergyLoss.hh.

Referenced by SetRndmStep().

G4double G4hRDEnergyLoss::RTable [static, protected]

Definition at line 181 of file G4hRDEnergyLoss.hh.

G4PhysicsTable * G4hRDEnergyLoss::theDEDXpbarTable = 0 [static, protected]

Definition at line 145 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable().

G4PhysicsTable * G4hRDEnergyLoss::theDEDXpTable = 0 [static, protected]

Definition at line 144 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable(), and G4hImpactIonisation::BuildPhysicsTable().

G4PhysicsTable * G4hRDEnergyLoss::theInverseRangepbarTable = 0 [static, protected]

Definition at line 151 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable().

G4PhysicsTable * G4hRDEnergyLoss::theInverseRangepTable = 0 [static, protected]

Definition at line 150 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable(), and G4hImpactIonisation::BuildPhysicsTable().

G4PhysicsTable * G4hRDEnergyLoss::theLabTimepbarTable = 0 [static, protected]

Definition at line 155 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable().

G4PhysicsTable * G4hRDEnergyLoss::theLabTimepTable = 0 [static, protected]

Definition at line 154 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable(), and G4hImpactIonisation::BuildPhysicsTable().

G4PhysicsTable* G4hRDEnergyLoss::theLossTable [protected]

Definition at line 185 of file G4hRDEnergyLoss.hh.

Referenced by G4hImpactIonisation::BuildPhysicsTable(), and ~G4hRDEnergyLoss().

G4PhysicsTable * G4hRDEnergyLoss::theProperTimepbarTable = 0 [static, protected]

Definition at line 158 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable().

G4PhysicsTable * G4hRDEnergyLoss::theProperTimepTable = 0 [static, protected]

Definition at line 157 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable(), and G4hImpactIonisation::BuildPhysicsTable().

G4PhysicsTable * G4hRDEnergyLoss::theRangepbarTable = 0 [static, protected]

Definition at line 147 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable().

G4PhysicsTable * G4hRDEnergyLoss::theRangepTable = 0 [static, protected]

Definition at line 146 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable(), and G4hImpactIonisation::BuildPhysicsTable().

G4int G4hRDEnergyLoss::TotBin = 360 [static, protected]

Definition at line 178 of file G4hRDEnergyLoss.hh.

Referenced by BuildDEDXTable(), and G4hImpactIonisation::BuildPhysicsTable().


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